CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/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(CLHEP::HepSymMatrix(6,0)); //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   static TrajectoryStateTransform transformer;
00060   //may be you want to do more than that if track does not go to IP
00061   return transformer. initialFreeState(track,thePropagator->magneticField());
00062 }
00063 
00064 double TrackAssociatorByPosition::quality(const TrajectoryStateOnSurface & tr, const TrajectoryStateOnSurface & sim) const {
00065   switch(theMethod){
00066   case 0:
00067     {
00068       AlgebraicVector5 v(tr.localParameters().vector() - sim.localParameters().vector());
00069       AlgebraicSymMatrix55 m(tr.localError().matrix());
00070       int ierr = ! m.Invert();
00071       if (ierr!=0) edm::LogInfo("TrackAssociatorByPosition")<<"error inverting the error matrix:\n"<<m;
00072       double est = ROOT::Math::Similarity(v,m);
00073       return est;
00074       break;
00075     }
00076   case 1:
00077     {
00078       return (tr.globalPosition() - sim.globalPosition()).mag();
00079       break;
00080     }
00081   case 2:
00082     {
00083       return  (deltaR<double>(tr.globalDirection().eta(),tr.globalDirection().phi(),
00084                               sim.globalDirection().eta(),sim.globalDirection().phi()));
00085       break;
00086     }
00087   case 3:
00088     {
00089       return (deltaR<double>(tr.globalPosition().eta(),tr.globalPosition().phi(),
00090                              sim.globalPosition().eta(),sim.globalPosition().phi()));
00091       break;
00092     }
00093   }
00094   //should never be reached
00095   edm::LogError("TrackAssociatorByPosition")<<"option: "<<theMethod<<" has not been recognized. association has no meaning.";
00096   return -1;
00097 }
00098 
00099 
00100 RecoToSimCollection TrackAssociatorByPosition::associateRecoToSim(const edm::RefToBaseVector<reco::Track>& tCH, 
00101                                                                   const edm::RefVector<TrackingParticleCollection>& tPCH,
00102                                                                   const edm::Event * e,
00103                                                                   const edm::EventSetup *setup ) const{
00104   RecoToSimCollection  outputCollection;
00105   //for each reco track find a matching tracking particle
00106   std::pair<unsigned int,unsigned int> minPair;
00107   const double dQmin_default=1542543;
00108   double dQmin=dQmin_default;
00109   for (unsigned int Ti=0; Ti!=tCH.size();++Ti){
00110     //initial state (initial OR inner OR outter)
00111     FreeTrajectoryState iState = getState(*(tCH)[Ti]);
00112 
00113     bool atLeastOne=false;
00114     //    for each tracking particle, find a state position and the plane to propagate the track to.
00115     for (unsigned int TPi=0;TPi!=tPCH.size();++TPi) {
00116       //get a state in the muon system 
00117       TrajectoryStateOnSurface simReferenceState = getState(*(tPCH)[TPi]);
00118       if (!simReferenceState.isValid()) continue;
00119 
00120       //propagate the TRACK to the surface
00121       TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface());
00122       if (!trackReferenceState.isValid()) continue; 
00123       
00124       //comparison
00125       double dQ= quality(trackReferenceState,simReferenceState);
00126       if (dQ < theQCut){
00127         atLeastOne=true;
00128         outputCollection.insert(tCH[Ti],
00129                                 std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH,TPi),-dQ));//association map with quality, is order greater-first
00130         edm::LogVerbatim("TrackAssociatorByPosition")<<"track number: "<<Ti
00131                                                      <<" associated with dQ: "<<dQ
00132                                                      <<" to TrackingParticle number: " <<TPi;}
00133       if (dQ < dQmin){
00134         dQmin=dQ;
00135         minPair = std::make_pair(Ti,TPi);}
00136     }//loop over tracking particles
00137     if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){
00138       outputCollection.insert(tCH[minPair.first],
00139                               std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH,minPair.second),-dQmin));}
00140   }//loop over tracks
00141   outputCollection.post_insert();
00142   return outputCollection;
00143 }
00144 
00145 
00146 
00147 SimToRecoCollection TrackAssociatorByPosition::associateSimToReco(const edm::RefToBaseVector<reco::Track>& tCH, 
00148                                                                   const edm::RefVector<TrackingParticleCollection>& tPCH,
00149                                                                   const edm::Event * e,
00150                                                                   const edm::EventSetup *setup ) const {
00151   SimToRecoCollection  outputCollection;
00152   //for each tracking particle, find matching tracks.
00153 
00154   std::pair<unsigned int,unsigned int> minPair;
00155   const double dQmin_default=1542543;
00156   double dQmin=dQmin_default;
00157   for (unsigned int TPi=0;TPi!=tPCH.size();++TPi){
00158     //get a state in the muon system
00159     TrajectoryStateOnSurface simReferenceState= getState(*(tPCH)[TPi]);
00160       
00161     if (!simReferenceState.isValid()) continue; 
00162     bool atLeastOne=false;
00163     //  propagate every track from any state (initial, inner, outter) to the surface 
00164     //  and make the position test
00165     for (unsigned int Ti=0; Ti!=tCH.size();++Ti){
00166       //initial state
00167       FreeTrajectoryState iState = getState(*(tCH)[Ti]);
00168         
00169       //propagation to surface
00170       TrajectoryStateOnSurface trackReferenceState = thePropagator->propagate(iState,simReferenceState.surface());
00171       if (!trackReferenceState.isValid()) continue;
00172         
00173       //comparison
00174       double dQ= quality(trackReferenceState, simReferenceState);
00175       if (dQ < theQCut){
00176         atLeastOne=true;
00177         outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH,TPi),
00178                                 std::make_pair(tCH[Ti],-dQ));//association map with quality, is order greater-first
00179         edm::LogVerbatim("TrackAssociatorByPosition")<<"TrackingParticle number: "<<TPi
00180                                                      <<" associated with dQ: "<<dQ
00181                                                      <<" to track number: "<<Ti;}
00182       if (dQ <dQmin){
00183         dQmin=dQ;
00184         minPair = std::make_pair(TPi,Ti);}
00185     }//loop over tracks
00186     if (theMinIfNoMatch && !atLeastOne && dQmin!=dQmin_default){
00187       outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH,minPair.first),
00188                               std::make_pair(tCH[minPair.second],-dQmin));}
00189   }//loop over tracking particles
00190   
00191   outputCollection.post_insert();
00192   return outputCollection;
00193 }