CMS 3D CMS Logo

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