CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/DataFormats/MuonReco/src/MuonCocktails.cc

Go to the documentation of this file.
00001 #include "TMath.h"
00002 #include "DataFormats/MuonReco/interface/MuonCocktails.h"
00003 #include "DataFormats/TrackReco/interface/Track.h"
00004 
00005 //
00006 // Return the TeV-optimized refit track (aka the cocktail or Tune P) or
00007 // the tracker track if either the optimized pT or tracker pT is below the pT threshold
00008 //
00009 reco::Muon::MuonTrackTypePair  muon::tevOptimized(const reco::TrackRef& combinedTrack,
00010                                                   const reco::TrackRef& trackerTrack,
00011                                                   const reco::TrackRef& tpfmsTrack,
00012                                                   const reco::TrackRef& pickyTrack,
00013                                                   const double ptThreshold,
00014                                                   const double tune1,
00015                                                   const double tune2) {
00016   
00017   // If Tracker pT is below the pT threshold (currently 200 GeV) - return the Tracker track
00018   if (trackerTrack->pt() < ptThreshold) return make_pair(trackerTrack,reco::Muon::InnerTrack);  
00019   
00020   // Array for convenience below.
00021   const reco::Muon::MuonTrackTypePair refit[4] = { 
00022     make_pair(trackerTrack, reco::Muon::InnerTrack), 
00023     make_pair(combinedTrack,reco::Muon::CombinedTrack),
00024     make_pair(tpfmsTrack,   reco::Muon::TPFMS),
00025     make_pair(pickyTrack,   reco::Muon::Picky)
00026   }; 
00027   
00028   // Calculate the log(tail probabilities). If there's a problem,
00029   // signify this with prob == 0. The current problems recognized are:
00030   // the track being not available, whether the (re)fit failed or it's
00031   // just not in the event, or if the (re)fit ended up with no valid
00032   // hits.
00033   double prob[4] = {0.};
00034   for (unsigned int i = 0; i < 4; ++i) 
00035     if (refit[i].first.isNonnull() && refit[i].first->numberOfValidHits()) 
00036       prob[i] = muon::trackProbability(refit[i].first); 
00037 
00038   //std::cout << "Probabilities: " << prob[0] << " " << prob[1] << " " << prob[2] << " " << prob[3] << std::endl;
00039   
00040   // Start with picky.
00041   int chosen = 3;
00042   
00043   // If there's a problem with picky, make the default one of the
00044   // other tracks. Try TPFMS first, then global, then tracker-only.
00045   if (prob[3] == 0.) { 
00046     if      (prob[2] > 0.) chosen = 2;
00047     else if (prob[1] > 0.) chosen = 1;
00048     else if (prob[0] > 0.) chosen = 0;
00049   } 
00050   
00051   // Now the algorithm: switch from picky to tracker-only if the
00052   // difference, log(tail prob(picky)) - log(tail prob(tracker-only))
00053   // is greater than a tuned value. Then compare the
00054   // so-picked track to TPFMS in the same manner using another tuned
00055   // value.
00056   if (prob[0] > 0. && prob[3] > 0. && (prob[3] - prob[0]) > tune1)
00057     chosen = 0;
00058   if (prob[2] > 0. && (prob[chosen] - prob[2]) > tune2)
00059     chosen = 2;
00060 
00061   // Done. If pT of the chosen track is below the threshold value, return the tracker track.
00062   if (refit[chosen].first->pt() < ptThreshold) return make_pair(trackerTrack,reco::Muon::InnerTrack);    
00063   
00064   // Return the chosen track (which can be the global track in
00065   // very rare cases).
00066   return refit[chosen];
00067 }
00068 
00069 //
00070 // calculate the tail probability (-ln(P)) of a fit
00071 //
00072 double muon::trackProbability(const reco::TrackRef track) {
00073   
00074   int nDOF = (int)track->ndof();
00075   if ( nDOF > 0 && track->chi2()> 0) { 
00076     return -log(TMath::Prob(track->chi2(), nDOF));
00077   } else { 
00078     return 0.0;
00079   }
00080   
00081 }
00082 
00083 reco::TrackRef muon::getTevRefitTrack(const reco::TrackRef& combinedTrack,
00084                                                      const reco::TrackToTrackMap& map) {
00085   reco::TrackToTrackMap::const_iterator it = map.find(combinedTrack);
00086   return it == map.end() ? reco::TrackRef() : it->val;
00087 }
00088 
00089 
00090 //
00091 // Get the sigma-switch decision (tracker-only versus global).
00092 //
00093 reco::Muon::MuonTrackTypePair muon::sigmaSwitch(const reco::TrackRef& combinedTrack,
00094                                                 const reco::TrackRef& trackerTrack,
00095                                                 const double nSigma,
00096                                                 const double ptThreshold) {
00097   // If either the global or tracker-only fits have pT below threshold
00098   // (default 200 GeV), return the tracker-only fit.
00099   if (combinedTrack->pt() < ptThreshold || trackerTrack->pt() < ptThreshold)
00100     return make_pair(trackerTrack,reco::Muon::InnerTrack);
00101   
00102   // If both are above the pT threshold, compare the difference in
00103   // q/p: if less than two sigma of the tracker-only track, switch to
00104   // global. Otherwise, use tracker-only.
00105   const double delta = fabs(trackerTrack->qoverp() - combinedTrack->qoverp());
00106   const double threshold = nSigma * trackerTrack->qoverpError();
00107   return delta > threshold ? make_pair(trackerTrack,reco::Muon::InnerTrack) :  make_pair(combinedTrack,reco::Muon::CombinedTrack);
00108 }
00109 
00110 //
00111 // Get the TMR decision (tracker-only versus TPFMS).
00112 //
00113 reco::Muon::MuonTrackTypePair muon::TMR(const reco::TrackRef& trackerTrack,
00114                                         const reco::TrackRef& fmsTrack,
00115                                         const double tune) {
00116   double probTK  = 0;
00117   double probFMS = 0;
00118   
00119   if (trackerTrack.isNonnull() && trackerTrack->numberOfValidHits())  
00120     probTK = muon::trackProbability(trackerTrack);
00121   if (fmsTrack.isNonnull() && fmsTrack->numberOfValidHits())
00122     probFMS = muon::trackProbability(fmsTrack);
00123   
00124   bool TKok  = probTK > 0;
00125   bool FMSok = probFMS > 0;
00126   
00127   if (TKok && FMSok) {
00128     if (probFMS - probTK > tune)
00129       return make_pair(trackerTrack,reco::Muon::InnerTrack);
00130     else
00131       return  make_pair(fmsTrack,reco::Muon::TPFMS);
00132   }
00133   else if (FMSok)
00134     return  make_pair(fmsTrack,reco::Muon::TPFMS);
00135   else if (TKok)
00136     return make_pair(trackerTrack,reco::Muon::InnerTrack);
00137   else
00138     return make_pair(reco::TrackRef(),reco::Muon::None); 
00139 }