CMS 3D CMS Logo

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