CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoMuon/GlobalTrackingTools/plugins/GlobalTrackQualityProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    GlobalTrackingTools
00004 // Class:      GlobalTrackQualityProducer
00005 //
00006 //
00007 // Original Author:  Adam Everett
00008 // $Id: GlobalTrackQualityProducer.cc,v 1.9 2013/01/06 19:16:52 dlange Exp $
00009 //
00010 //
00011 
00012 // system include files
00013 #include <memory>
00014 
00015 // user include files
00016 #include "FWCore/Framework/interface/Frameworkfwd.h"
00017 #include "FWCore/Framework/interface/EDProducer.h"
00018 #include "FWCore/Utilities/interface/isFinite.h"
00019 
00020 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00021 #include "RecoMuon/GlobalTrackingTools/plugins/GlobalTrackQualityProducer.h"
00022 
00023 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00024 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00025 #include "DataFormats/MuonReco/interface/MuonQuality.h"
00026 #include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
00027 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00028 
00029 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00030 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00031 
00032 GlobalTrackQualityProducer::GlobalTrackQualityProducer(const edm::ParameterSet& iConfig):
00033   inputCollection_(iConfig.getParameter<edm::InputTag>("InputCollection")),inputLinksCollection_(iConfig.getParameter<edm::InputTag>("InputLinksCollection")),theService(0),theGlbRefitter(0),theGlbMatcher(0)
00034 {
00035   // service parameters
00036   edm::ParameterSet serviceParameters = iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
00037   theService = new MuonServiceProxy(serviceParameters);     
00038   
00039   // TrackRefitter parameters
00040   edm::ParameterSet refitterParameters = iConfig.getParameter<edm::ParameterSet>("RefitterParameters");
00041   theGlbRefitter = new GlobalMuonRefitter(refitterParameters, theService);
00042 
00043   edm::ParameterSet trackMatcherPSet = iConfig.getParameter<edm::ParameterSet>("GlobalMuonTrackMatcher");
00044   theGlbMatcher = new GlobalMuonTrackMatcher(trackMatcherPSet,theService);
00045 
00046   double maxChi2 = iConfig.getParameter<double>("MaxChi2");
00047   double nSigma = iConfig.getParameter<double>("nSigma");
00048   theEstimator = new Chi2MeasurementEstimator(maxChi2,nSigma);
00049  
00050   produces<edm::ValueMap<reco::MuonQuality> >();
00051 }
00052 
00053 GlobalTrackQualityProducer::~GlobalTrackQualityProducer() {
00054   if (theService) delete theService;
00055   if (theGlbRefitter) delete theGlbRefitter;
00056   if (theGlbMatcher) delete theGlbMatcher;
00057   if (theEstimator) delete theEstimator;
00058 }
00059 
00060 void
00061 GlobalTrackQualityProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00062 {
00063   const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
00064   
00065   theService->update(iSetup);
00066 
00067   theGlbRefitter->setEvent(iEvent);
00068 
00069   theGlbRefitter->setServices(theService->eventSetup());
00070 
00071   // Take the GLB muon container(s)
00072   edm::Handle<reco::TrackCollection> glbMuons;
00073   iEvent.getByLabel(inputCollection_,glbMuons);
00074   
00075   edm::Handle<reco::MuonTrackLinksCollection>    linkCollectionHandle;
00076   iEvent.getByLabel(inputLinksCollection_, linkCollectionHandle);
00077 
00078   //Retrieve tracker topology from geometry
00079   edm::ESHandle<TrackerTopology> tTopoHand;
00080   iSetup.get<IdealGeometryRecord>().get(tTopoHand);
00081   const TrackerTopology *tTopo=tTopoHand.product();
00082 
00083 
00084   // reserve some space
00085   std::vector<reco::MuonQuality> valuesQual;
00086   valuesQual.reserve(glbMuons->size());
00087   
00088   int trackIndex = 0;
00089   for (reco::TrackCollection::const_iterator track = glbMuons->begin(); track!=glbMuons->end(); ++track , ++trackIndex) {
00090     reco::TrackRef glbRef(glbMuons,trackIndex);
00091     reco::TrackRef staTrack = reco::TrackRef();
00092 
00093     std::vector<Trajectory> refitted=theGlbRefitter->refit(*track,1,tTopo);
00094 
00095     LogTrace(theCategory)<<"GLBQual N refitted " << refitted.size();
00096     
00097     std::pair<double,double> thisKink;
00098     double relative_muon_chi2 = 0.0;
00099     double relative_tracker_chi2 = 0.0;
00100     double glbTrackProbability = 0.0;
00101     if(refitted.size()>0) {
00102       thisKink = kink(refitted.front()) ;      
00103       std::pair<double,double> chi = newChi2(refitted.front());
00104       relative_muon_chi2 = chi.second; //normalized inside to /sum(muHits.dimension)
00105       relative_tracker_chi2 = chi.first; //normalized inside to /sum(tkHits.dimension)
00106       glbTrackProbability = trackProbability(refitted.front());
00107     }
00108 
00109     LogTrace(theCategory)<<"GLBQual: Kink " << thisKink.first << " " << thisKink.second;
00110     LogTrace(theCategory)<<"GLBQual: Rel Chi2 " << relative_tracker_chi2 << " " << relative_muon_chi2;
00111     LogTrace(theCategory)<<"GLBQual: trackProbability " << glbTrackProbability;
00112 
00113     // Fill the STA-TK match information
00114     float chi2, d, dist, Rpos;
00115     chi2 = d = dist = Rpos = -1.0;
00116     bool passTight = false;
00117     typedef MuonTrajectoryBuilder::TrackCand TrackCand;
00118     if ( linkCollectionHandle.isValid() ) {
00119       for ( reco::MuonTrackLinksCollection::const_iterator links = linkCollectionHandle->begin();
00120             links != linkCollectionHandle->end(); ++links )
00121         {
00122           if ( links->trackerTrack().isNull() ||
00123                links->standAloneTrack().isNull() ||
00124                links->globalTrack().isNull() ) 
00125             {
00126               edm::LogWarning(theCategory) << "Global muon links to constituent tracks are invalid. There should be no such object. Muon is skipped.";
00127               continue;
00128             }
00129           if (links->globalTrack() == glbRef) {
00130             staTrack = !links->standAloneTrack().isNull() ? links->standAloneTrack() : reco::TrackRef();
00131             TrackCand staCand = TrackCand((Trajectory*)(0),links->standAloneTrack());
00132             TrackCand tkCand = TrackCand((Trajectory*)(0),links->trackerTrack());
00133             chi2 = theGlbMatcher->match(staCand,tkCand,0,0);
00134             d    = theGlbMatcher->match(staCand,tkCand,1,0);
00135             Rpos = theGlbMatcher->match(staCand,tkCand,2,0);
00136             dist = theGlbMatcher->match(staCand,tkCand,3,0);
00137             passTight = theGlbMatcher->matchTight(staCand,tkCand);
00138           }
00139         }
00140     }
00141 
00142     if(!staTrack.isNull()) LogTrace(theCategory)<<"GLBQual: Used UpdatedAtVtx : " <<  (iEvent.getProvenance(staTrack.id()).productInstanceName() == std::string("UpdatedAtVtx"));
00143 
00144     float maxFloat01 = std::numeric_limits<float>::max()*0.1; // a better solution would be to use float above .. m/be not
00145     reco::MuonQuality muQual;
00146     if(!staTrack.isNull()) muQual.updatedSta = iEvent.getProvenance(staTrack.id()).productInstanceName() == std::string("UpdatedAtVtx");
00147     muQual.trkKink    = thisKink.first > maxFloat01 ? maxFloat01 : thisKink.first;
00148     muQual.glbKink    = thisKink.second > maxFloat01 ? maxFloat01 : thisKink.second;
00149     muQual.trkRelChi2 = relative_tracker_chi2 > maxFloat01 ? maxFloat01 : relative_tracker_chi2;
00150     muQual.staRelChi2 = relative_muon_chi2 > maxFloat01 ? maxFloat01 : relative_muon_chi2;
00151     muQual.tightMatch = passTight;
00152     muQual.chi2LocalPosition = dist;
00153     muQual.chi2LocalMomentum = chi2;
00154     muQual.localDistance = d;
00155     muQual.globalDeltaEtaPhi = Rpos;
00156     muQual.glbTrackProbability = glbTrackProbability;
00157     valuesQual.push_back(muQual);
00158   }
00159 
00160   /*
00161   for(int i = 0; i < valuesTkRelChi2.size(); i++) {
00162     LogTrace(theCategory)<<"value " << valuesTkRelChi2[i] ;
00163   }
00164   */
00165 
00166   // create and fill value maps
00167   std::auto_ptr<edm::ValueMap<reco::MuonQuality> > outQual(new edm::ValueMap<reco::MuonQuality>());
00168   edm::ValueMap<reco::MuonQuality>::Filler fillerQual(*outQual);
00169   fillerQual.insert(glbMuons, valuesQual.begin(), valuesQual.end());
00170   fillerQual.fill();
00171   
00172   // put value map into event
00173   iEvent.put(outQual);
00174 }
00175 
00176 std::pair<double,double> GlobalTrackQualityProducer::kink(Trajectory& muon) const {
00177 
00178   const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
00179    
00180   using namespace std;
00181   using namespace edm;
00182   using namespace reco;
00183  
00184   double result = 0.0;
00185   double resultGlb = 0.0;
00186 
00187      
00188   typedef TransientTrackingRecHit::ConstRecHitPointer   ConstRecHitPointer;
00189   typedef ConstRecHitPointer RecHit;
00190   typedef std::vector<TrajectoryMeasurement>::const_iterator TMI;
00191 
00192 
00193   vector<TrajectoryMeasurement> meas = muon.measurements();
00194   
00195   for ( TMI m = meas.begin(); m != meas.end(); m++ ) {
00196     TransientTrackingRecHit::ConstRecHitPointer hit = m->recHit();
00197 
00198     //not used    double estimate = 0.0;
00199 
00200     RecHit rhit = (*m).recHit();
00201     bool ok = false;
00202     if ( rhit->isValid() ) {
00203       if(DetId::Tracker == rhit->geographicalId().det()) ok = true;
00204     }
00205 
00206     //if ( !ok ) continue;
00207     
00208     const TrajectoryStateOnSurface& tsos = (*m).predictedState();
00209 
00210 
00211     if ( tsos.isValid() && rhit->isValid() && rhit->hit()->isValid()
00212          && !edm::isNotFinite(rhit->localPositionError().xx()) //this is paranoia induced by reported case
00213          && !edm::isNotFinite(rhit->localPositionError().xy()) //it's better to track down the origin of bad numbers
00214          && !edm::isNotFinite(rhit->localPositionError().yy())
00215          ) {
00216 
00217       double phi1 = tsos.globalPosition().phi();
00218       if ( phi1 < 0 ) phi1 = 2*M_PI + phi1;
00219 
00220       double phi2 = rhit->globalPosition().phi();
00221       if ( phi2 < 0 ) phi2 = 2*M_PI + phi2;
00222 
00223       double diff = fabs(phi1 - phi2);
00224       if ( diff > M_PI ) diff = 2*M_PI - diff;
00225 
00226       GlobalPoint hitPos = rhit->globalPosition();
00227 
00228       GlobalError hitErr = rhit->globalPositionError();
00229       //LogDebug(theCategory)<<"hitPos " << hitPos;
00230       double error = hitErr.phierr(hitPos);  // error squared
00231 
00232       double s = ( error > 0.0 ) ? (diff*diff)/error : (diff*diff);
00233 
00234       if(ok) result += s;
00235       resultGlb += s;
00236     }
00237     
00238   }
00239   
00240   
00241   return std::pair<double,double>(result,resultGlb);
00242   
00243 }
00244 
00245 std::pair<double,double> GlobalTrackQualityProducer::newChi2(Trajectory& muon) const {
00246   const std::string theCategory = "Muon|RecoMuon|GlobalTrackQualityProducer";
00247 
00248   using namespace std;
00249   using namespace edm;
00250   using namespace reco;
00251 
00252   double muChi2 = 0.0;
00253   double tkChi2 = 0.0;
00254   unsigned int muNdof = 0;
00255   unsigned int tkNdof = 0;
00256   
00257   typedef TransientTrackingRecHit::ConstRecHitPointer   ConstRecHitPointer;
00258   typedef ConstRecHitPointer RecHit;
00259   typedef vector<TrajectoryMeasurement>::const_iterator TMI;
00260 
00261   vector<TrajectoryMeasurement> meas = muon.measurements();
00262 
00263   for ( TMI m = meas.begin(); m != meas.end(); m++ ) {
00264     TransientTrackingRecHit::ConstRecHitPointer hit = m->recHit();
00265     const TrajectoryStateOnSurface& uptsos = (*m).updatedState();
00266     TransientTrackingRecHit::RecHitPointer preciseHit = hit->clone(uptsos);
00267     double estimate = 0.0;
00268     if (preciseHit->isValid() && uptsos.isValid()) {
00269       estimate = theEstimator->estimate(uptsos, *preciseHit ).second;
00270     }
00271     
00272     //LogTrace(theCategory) << "estimate " << estimate << " TM.est " << m->estimate();
00273     //UNUSED:    double tkDiff = 0.0;
00274     //UNUSED:    double staDiff = 0.0;
00275     if ( hit->isValid() &&  (hit->geographicalId().det()) == DetId::Tracker ) {
00276       tkChi2 += estimate;
00277       //UNUSED:      tkDiff = estimate - m->estimate();
00278       tkNdof += hit->dimension();
00279     }
00280     if ( hit->isValid() &&  (hit->geographicalId().det()) == DetId::Muon ) {
00281       muChi2 += estimate;
00282       //UNUSED      staDiff = estimate - m->estimate();
00283       muNdof += hit->dimension();
00284     }
00285   }
00286   
00287   if (tkNdof < 6 ) tkChi2 = tkChi2; // or should I set it to  a large number ?
00288   else tkChi2 /= (tkNdof-5.);
00289 
00290   if (muNdof < 6 ) muChi2 = muChi2; // or should I set it to  a large number ?
00291   else muChi2 /= (muNdof-5.);
00292 
00293   return std::pair<double,double>(tkChi2,muChi2);
00294        
00295 }
00296 
00297 //
00298 // calculate the tail probability (-ln(P)) of a fit
00299 //
00300 double 
00301 GlobalTrackQualityProducer::trackProbability(Trajectory& track) const {
00302 
00303   if ( track.ndof() > 0 && track.chiSquared() > 0 ) { 
00304     return -LnChiSquaredProbability(track.chiSquared(), track.ndof());
00305   } else {
00306     return 0.0;
00307   }
00308 
00309 }
00310 
00311 //#include "FWCore/Framework/interface/MakerMacros.h"
00312 //DEFINE_FWK_MODULE(GlobalTrackQualityProducer);