CMS 3D CMS Logo

Classes | Enumerations | Functions

muon Namespace Reference

Classes

struct  SelectionTypeStringToEnum
 a lightweight "map" for selection type string label and enum value More...

Enumerations

enum  AlgorithmType { TMLastStation, TM2DCompatibility, TMOneStation }
enum  SelectionType {
  All = 0, AllGlobalMuons = 1, AllStandAloneMuons = 2, AllTrackerMuons = 3,
  TrackerMuonArbitrated = 4, AllArbitrated = 5, GlobalMuonPromptTight = 6, TMLastStationLoose = 7,
  TMLastStationTight = 8, TM2DCompatibilityLoose = 9, TM2DCompatibilityTight = 10, TMOneStationLoose = 11,
  TMOneStationTight = 12, TMLastStationOptimizedLowPtLoose = 13, TMLastStationOptimizedLowPtTight = 14, GMTkChiCompatibility = 15,
  GMStaChiCompatibility = 16, GMTkKinkTight = 17, TMLastStationAngLoose = 18, TMLastStationAngTight = 19,
  TMOneStationAngLoose = 20, TMOneStationAngTight = 21, TMLastStationOptimizedBarrelLowPtLoose = 22, TMLastStationOptimizedBarrelLowPtTight = 23
}
 

Selector type.

More...

Functions

float caloCompatibility (const reco::Muon &muon)
reco::TrackRef getTevRefitTrack (const reco::TrackRef &combinedTrack, const reco::TrackToTrackMap &map)
bool isGoodMuon (const reco::Muon &muon, SelectionType type, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
 main GoodMuon wrapper call
bool isGoodMuon (const reco::Muon &muon, AlgorithmType type, double minCompatibility, reco::Muon::ArbitrationType arbitrationType)
bool isGoodMuon (const reco::Muon &muon, AlgorithmType type, int minNumberOfMatches, double maxAbsDx, double maxAbsPullX, double maxAbsDy, double maxAbsPullY, double maxChamberDist, double maxChamberDistPull, reco::Muon::ArbitrationType arbitrationType, bool syncMinNMatchesNRequiredStationsInBarrelOnly=true, bool applyAlsoAngularCuts=false)
bool isTightMuon (const reco::Muon &, const reco::Vertex &)
bool overlap (const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
unsigned int RequiredStationMask (const reco::Muon &muon, double maxChamberDist, double maxChamberDistPull, reco::Muon::ArbitrationType arbitrationType)
float segmentCompatibility (const reco::Muon &muon, reco::Muon::ArbitrationType arbitrationType=reco::Muon::SegmentAndTrackArbitration)
SelectionType selectionTypeFromString (const std::string &label)
reco::Muon::MuonTrackTypePair sigmaSwitch (const reco::TrackRef &combinedTrack, const reco::TrackRef &trackerTrack, const double nSigma=2., const double ptThreshold=200.)
reco::Muon::MuonTrackTypePair sigmaSwitch (const reco::Muon &muon, const double nSigma=2., const double ptThreshold=200.)
reco::Muon::MuonTrackTypePair tevOptimized (const reco::Muon &muon, const reco::TrackToTrackMap &tevMap1, const reco::TrackToTrackMap &tevMap2, const reco::TrackToTrackMap &tevMap3)
reco::Muon::MuonTrackTypePair tevOptimized (const reco::TrackRef &combinedTrack, const reco::TrackRef &trackerTrack, const reco::TrackToTrackMap &tevMap1, const reco::TrackToTrackMap &tevMap2, const reco::TrackToTrackMap &tevMap3, const double ptThreshold=200., const double tune1=4., const double tune2=6.)
reco::Muon::MuonTrackTypePair tevOptimized (const reco::Muon &muon, const double ptThreshold=200., const double tune1=4., const double tune2=6.)
reco::Muon::MuonTrackTypePair tevOptimized (const reco::TrackRef &combinedTrack, const reco::TrackRef &trackerTrack, const reco::TrackRef &tpfmsTrack, const reco::TrackRef &pickyTrack, const double ptThreshold=200., const double tune1=4., const double tune2=6.)
reco::Muon::MuonTrackTypePair TMR (const reco::TrackRef &trackerTrack, const reco::TrackRef &fmsTrack, const double tune=4.)
double trackProbability (const reco::TrackRef track)

Enumeration Type Documentation

Enumerator:
TMLastStation 
TM2DCompatibility 
TMOneStation 

Definition at line 63 of file MuonSelectors.h.

Selector type.

Enumerator:
All 
AllGlobalMuons 
AllStandAloneMuons 
AllTrackerMuons 
TrackerMuonArbitrated 
AllArbitrated 
GlobalMuonPromptTight 
TMLastStationLoose 
TMLastStationTight 
TM2DCompatibilityLoose 
TM2DCompatibilityTight 
TMOneStationLoose 
TMOneStationTight 
TMLastStationOptimizedLowPtLoose 
TMLastStationOptimizedLowPtTight 
GMTkChiCompatibility 
GMStaChiCompatibility 
GMTkKinkTight 
TMLastStationAngLoose 
TMLastStationAngTight 
TMOneStationAngLoose 
TMOneStationAngTight 
TMLastStationOptimizedBarrelLowPtLoose 
TMLastStationOptimizedBarrelLowPtTight 

Definition at line 18 of file MuonSelectors.h.

                      {
      All = 0,                      // dummy options - always true
      AllGlobalMuons = 1,           // checks isGlobalMuon flag
      AllStandAloneMuons = 2,       // checks isStandAloneMuon flag
      AllTrackerMuons = 3,          // checks isTrackerMuon flag
      TrackerMuonArbitrated = 4,    // resolve ambiguity of sharing segments
      AllArbitrated = 5,            // all muons with the tracker muon arbitrated
      GlobalMuonPromptTight = 6,    // global muons with tighter fit requirements
      TMLastStationLoose = 7,       // penetration depth loose selector
      TMLastStationTight = 8,       // penetration depth tight selector
      TM2DCompatibilityLoose = 9,   // likelihood based loose selector
      TM2DCompatibilityTight = 10,  // likelihood based tight selector
      TMOneStationLoose = 11,       // require one well matched segment
      TMOneStationTight = 12,       // require one well matched segment
      TMLastStationOptimizedLowPtLoose = 13, // combination of TMLastStation and TMOneStation
      TMLastStationOptimizedLowPtTight = 14, // combination of TMLastStation and TMOneStation
      GMTkChiCompatibility = 15,    // require tk stub have good chi2 relative to glb track
      GMStaChiCompatibility = 16,   // require sta stub have good chi2 compatibility relative to glb track
      GMTkKinkTight = 17,           // require a small kink value in the tracker stub
      TMLastStationAngLoose = 18,   // TMLastStationLoose with additional angular cuts
      TMLastStationAngTight = 19,   // TMLastStationTight with additional angular cuts
      TMOneStationAngLoose = 20,    // TMOneStationLoose with additional angular cuts
      TMOneStationAngTight = 21,    // TMOneStationTight with additional angular cuts
      // The two algorithms that follow are identical to what were known as
      // TMLastStationOptimizedLowPt* (sans the Barrel) as late as revision
      // 1.7 of this file. The names were changed because indeed the low pt
      // optimization applies only to the barrel region, whereas the sel-
      // ectors above are more efficient at low pt in the endcaps, which is
      // what we feel is more suggestive of the algorithm name. This will be
      // less confusing for future generations of CMS members, I hope...
      TMLastStationOptimizedBarrelLowPtLoose = 22, // combination of TMLastStation and TMOneStation but with low pT optimization in barrel only
      TMLastStationOptimizedBarrelLowPtTight = 23  // combination of TMLastStation and TMOneStation but with low pT optimization in barrel only
   };

Function Documentation

float muon::caloCompatibility ( const reco::Muon muon)
reco::TrackRef muon::getTevRefitTrack ( const reco::TrackRef combinedTrack,
const reco::TrackToTrackMap map 
)

Definition at line 83 of file MuonCocktails.cc.

References edm::AssociationMap< Tag >::end(), and edm::AssociationMap< Tag >::find().

Referenced by MuonIdProducer::makeMuon(), and tevOptimized().

                                                                                     {
  reco::TrackToTrackMap::const_iterator it = map.find(combinedTrack);
  return it == map.end() ? reco::TrackRef() : it->val;
}
bool muon::isGoodMuon ( const reco::Muon muon,
SelectionType  type,
reco::Muon::ArbitrationType  arbitrationType = reco::Muon::SegmentAndTrackArbitration 
)

main GoodMuon wrapper call

Definition at line 534 of file MuonSelectors.cc.

References All, AllArbitrated, AllGlobalMuons, AllStandAloneMuons, AllTrackerMuons, reco::Muon::combinedQuality(), reco::LeafCandidate::eta(), GlobalMuonPromptTight, reco::Muon::globalTrack(), GMStaChiCompatibility, GMTkChiCompatibility, GMTkKinkTight, reco::Muon::innerTrack(), reco::Muon::isGlobalMuon(), reco::Muon::isQualityValid(), reco::Muon::isStandAloneMuon(), reco::Muon::isTrackerMuon(), reco::Muon::numberOfMatches(), reco::Muon::outerTrack(), reco::LeafCandidate::pt(), TM2DCompatibility, TM2DCompatibilityLoose, TM2DCompatibilityTight, TMLastStation, TMLastStationAngLoose, TMLastStationAngTight, TMLastStationLoose, TMLastStationOptimizedBarrelLowPtLoose, TMLastStationOptimizedBarrelLowPtTight, TMLastStationOptimizedLowPtLoose, TMLastStationOptimizedLowPtTight, TMLastStationTight, TMOneStation, TMOneStationAngLoose, TMOneStationAngTight, TMOneStationLoose, TMOneStationTight, and TrackerMuonArbitrated.

Referenced by ExampleMuonAnalyzer::analyze(), MuonIdVal::analyze(), TrackEfficiencyMonitor::analyze(), HiggsDQM::analyze(), BPhysicsOniaDQM::analyze(), MuonCosmicCompatibilityFiller::checkMuonID(), PFRecoTauDiscriminationAgainstMuon::discriminate(), WMuNuSelector::filter(), WMuNuValidator::filter(), HLTDiMuonGlbTrkFilter::hltFilter(), PFMuonAlgo::isGlobalLooseMuon(), PFMuonAlgo::isGlobalTightMuon(), isTightMuon(), PFMuonAlgo::isTightMuonPOG(), PFMuonAlgo::isTrackerLooseMuon(), PFMuonAlgo::isTrackerTightMuon(), JetPlusTrackCorrector::matchMuons(), pat::MuonSelector::muIdSelection_(), pat::Muon::muonID(), PFMuonAlgo::printMuonProperties(), SoftLepton::produce(), MuonTrackProducer::produce(), MuonRefProducer::produce(), MuonSelectionTypeValueMapProducer::produce(), and FourVectorHLTOffline::selectMuons().

{
  switch (type)
    {
    case muon::All:
      return true;
      break;
    case muon::AllGlobalMuons:
      return muon.isGlobalMuon();
      break;
    case muon::AllTrackerMuons:
      return muon.isTrackerMuon();
      break;
    case muon::AllStandAloneMuons:
      return muon.isStandAloneMuon();
      break;
    case muon::TrackerMuonArbitrated:
      return muon.isTrackerMuon() && muon.numberOfMatches(arbitrationType)>0;
      break;
    case muon::AllArbitrated:
      return ! muon.isTrackerMuon() || muon.numberOfMatches(arbitrationType)>0;
      break;
    case muon::GlobalMuonPromptTight:
      return muon.isGlobalMuon() && muon.globalTrack()->normalizedChi2()<10. && muon.globalTrack()->hitPattern().numberOfValidMuonHits() >0;
      break;
      // For "Loose" algorithms we choose maximum y quantity cuts of 1E9 instead of
      // 9999 as before.  We do this because the muon methods return 999999 (note
      // there are six 9's) when the requested information is not available.  For
      // example, if a muon fails to traverse the z measuring superlayer in a station
      // in the DT, then all methods involving segmentY in this station return
      // 999999 to demonstrate that the information is missing.  In order to not
      // penalize muons for missing y information in Loose algorithms where we do
      // not care at all about y information, we raise these limits.  In the
      // TMLastStation and TMOneStation algorithms we actually use this huge number
      // to determine whether to consider y information at all.
    case muon::TMLastStationLoose:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,true,false);
      break;
    case muon::TMLastStationTight:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,true,false);
      break;
    case muon::TMOneStationLoose:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
      break;
    case muon::TMOneStationTight:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
      break;
    case muon::TMLastStationOptimizedLowPtLoose:
      if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
        return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
      else
        return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,false,false);
      break;
    case muon::TMLastStationOptimizedLowPtTight:
      if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
        return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
      else
        return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,false,false);
      break;
      //compatibility loose
    case muon::TM2DCompatibilityLoose:
      return muon.isTrackerMuon() && isGoodMuon(muon,TM2DCompatibility,0.7,arbitrationType);
      break;
      //compatibility tight
    case muon::TM2DCompatibilityTight:
      return muon.isTrackerMuon() && isGoodMuon(muon,TM2DCompatibility,1.0,arbitrationType);
      break;
    case muon::GMTkChiCompatibility:
      return muon.isGlobalMuon() && muon.isQualityValid() && fabs(muon.combinedQuality().trkRelChi2 - muon.innerTrack()->normalizedChi2()) < 2.0;
      break;
    case muon::GMStaChiCompatibility:
      return muon.isGlobalMuon() && muon.isQualityValid() && fabs(muon.combinedQuality().staRelChi2 - muon.outerTrack()->normalizedChi2()) < 2.0;
      break;
    case muon::GMTkKinkTight:
      return muon.isGlobalMuon() && muon.isQualityValid() && muon.combinedQuality().trkKink < 100.0;
      break;
    case muon::TMLastStationAngLoose:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,false,true);
      break;
    case muon::TMLastStationAngTight:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,false,true);
      break;
    case muon::TMOneStationAngLoose:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,true);
      break;
    case muon::TMOneStationAngTight:
      return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,true);
      break;
    case muon::TMLastStationOptimizedBarrelLowPtLoose:
      if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
        return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
      else
        return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,true,false);
      break;
    case muon::TMLastStationOptimizedBarrelLowPtTight:
      if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
        return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
      else
        return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,true,false);
      break;
    default:
      return false;
    }
}
bool muon::isGoodMuon ( const reco::Muon muon,
AlgorithmType  type,
double  minCompatibility,
reco::Muon::ArbitrationType  arbitrationType 
)

Definition at line 287 of file MuonSelectors.cc.

References caloCompatibility(), reco::Muon::isMatchesValid(), segmentCompatibility(), and TM2DCompatibility.

                                                                   {
  if (!muon.isMatchesValid()) return false;
  bool goodMuon = false;
  
  switch( type ) {
  case TM2DCompatibility:
    // Simplistic first cut in the 2D segment- vs calo-compatibility plane. Will have to be refined!
    if( ( (0.8*caloCompatibility( muon ))+(1.2*segmentCompatibility( muon, arbitrationType )) ) > minCompatibility ) goodMuon = true;
    else goodMuon = false;
    return goodMuon;
    break;
  default : 
    //  LogTrace("MuonIdentification")<<"            // Invalid Algorithm Type called!";
    goodMuon = false;
    return goodMuon;
    break;
  }
}
bool muon::isGoodMuon ( const reco::Muon muon,
AlgorithmType  type,
int  minNumberOfMatches,
double  maxAbsDx,
double  maxAbsPullX,
double  maxAbsDy,
double  maxAbsPullY,
double  maxChamberDist,
double  maxChamberDistPull,
reco::Muon::ArbitrationType  arbitrationType,
bool  syncMinNMatchesNRequiredStationsInBarrelOnly = true,
bool  applyAlsoAngularCuts = false 
)

Definition at line 309 of file MuonSelectors.cc.

References reco::Muon::dX(), reco::Muon::dY(), reco::LeafCandidate::eta(), reco::Muon::isMatchesValid(), reco::Muon::pullDxDz(), reco::Muon::pullDyDz(), reco::Muon::pullX(), reco::Muon::pullY(), RequiredStationMask(), relativeConstraints::station, reco::Muon::stationMask(), TMLastStation, and TMOneStation.

{
   if (!muon.isMatchesValid()) return false;
   bool goodMuon = false;

   if (type == TMLastStation) {
      // To satisfy my own paranoia, if the user specifies that the
      // minimum number of matches is zero, then return true.
      if(minNumberOfMatches == 0) return true;

      unsigned int theStationMask = muon.stationMask(arbitrationType);
      unsigned int theRequiredStationMask = RequiredStationMask(muon, maxChamberDist, maxChamberDistPull, arbitrationType);

      // Require that there be at least a minimum number of segments
      int numSegs = 0;
      int numRequiredStations = 0;
      for(int it = 0; it < 8; ++it) {
         if(theStationMask & 1<<it) ++numSegs;
         if(theRequiredStationMask & 1<<it) ++numRequiredStations;
      }

      // Make sure the minimum number of matches is not greater than
      // the number of required stations but still greater than zero
      if (syncMinNMatchesNRequiredStationsInBarrelOnly) {
         // Note that we only do this in the barrel region!
         if (fabs(muon.eta()) < 1.2) {
            if(minNumberOfMatches > numRequiredStations)
               minNumberOfMatches = numRequiredStations;
            if(minNumberOfMatches < 1) //SK: this only happens for negative values
               minNumberOfMatches = 1;
         }
      } else {
         if(minNumberOfMatches > numRequiredStations)
            minNumberOfMatches = numRequiredStations;
         if(minNumberOfMatches < 1) //SK: this only happens for negative values
            minNumberOfMatches = 1;
      }

      if(numSegs >= minNumberOfMatches) goodMuon = 1;

      // Require that last required station have segment
      // If there are zero required stations keep track
      // of the last station with a segment so that we may
      // apply the quality cuts below to it instead
      int lastSegBit = 0;
      if(theRequiredStationMask) {
         for(int stationIdx = 7; stationIdx >= 0; --stationIdx)
           if(theRequiredStationMask & 1<<stationIdx){
               if(theStationMask & 1<<stationIdx) {
                  lastSegBit = stationIdx;
                  goodMuon &= 1;
                  break;
               } else {
                  goodMuon = false;
                  break;
               }
           }
      } else {
         for(int stationIdx = 7; stationIdx >= 0; --stationIdx)
            if(theStationMask & 1<<stationIdx) {
               lastSegBit = stationIdx;
               break;
            }
      }

      if(!goodMuon) return false;

      // Impose pull cuts on last segment
      int station = 0, detector = 0;
      station  = lastSegBit < 4 ? lastSegBit+1 : lastSegBit-3;
      detector = lastSegBit < 4 ? 1 : 2;

      // Check x information
      if(fabs(muon.pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
            fabs(muon.dX(station,detector,arbitrationType)) > maxAbsDx)
         return false;

      if(applyAlsoAngularCuts && fabs(muon.pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX)
         return false;

      // Is this a tight algorithm, i.e. do we bother to check y information?
      if (maxAbsDy < 999999) { // really if maxAbsDy < 1E9 as currently defined

         // Check y information
         if (detector == 2) { // CSC
            if(fabs(muon.pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
                  fabs(muon.dY(station,2,arbitrationType)) > maxAbsDy)
               return false;

            if(applyAlsoAngularCuts && fabs(muon.pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY)
               return false;
         } else {
            //
            // In DT, if this is a "Tight" algorithm and the last segment is
            // missing y information (always the case in station 4!!!), impose
            // respective cuts on the next station in the stationMask that has
            // a segment with y information.  If there are no segments with y
            // information then there is nothing to penalize. Should we
            // penalize in Tight for having zero segments with y information?
            // That is the fundamental question.  Of course I am being uber
            // paranoid; if this is a good muon then there will probably be at
            // least one segment with y information but not always.  Suppose
            // somehow a muon only creates segments in station 4, then we
            // definitely do not want to require that there be at least one
            // segment with y information because we will lose it completely.
            //

            for (int stationIdx = station; stationIdx > 0; --stationIdx) {
               if(! (theStationMask & 1<<(stationIdx-1)))  // don't bother if the station is not in the stationMask
                  continue;

               if(muon.dY(stationIdx,1,arbitrationType) > 999998) // no y-information
                  continue;

               if(fabs(muon.pullY(stationIdx,1,arbitrationType,1)) > maxAbsPullY &&
                     fabs(muon.dY(stationIdx,1,arbitrationType)) > maxAbsDy) {
                  return false;
               }

               if(applyAlsoAngularCuts && fabs(muon.pullDyDz(stationIdx,1,arbitrationType,1)) > maxAbsPullY)
                  return false;

               // If we get this far then great this is a good muon
               return true;
            }
         }
      }

      return goodMuon;
   } // TMLastStation

   // TMOneStation requires only that there be one "good" segment, regardless
   // of the required stations.  We do not penalize if there are absolutely zero
   // segments with y information in the Tight algorithm.  Maybe I'm being
   // paranoid but so be it.  If it's really a good muon then we will probably
   // find at least one segment with both x and y information but you never
   // know, and I don't want to deal with a potential inefficiency in the DT
   // like we did with the original TMLastStation.  Incidentally, not penalizing
   // for total lack of y information in the Tight algorithm is what is done in
   // the new TMLastStation
   //                   
   if (type == TMOneStation) {
      unsigned int theStationMask = muon.stationMask(arbitrationType);

      // Of course there must be at least one segment
      if (! theStationMask) return false;

      int  station = 0, detector = 0;
      // Keep track of whether or not there is a DT segment with y information.
      // In the end, if it turns out there are absolutely zero DT segments with
      // y information, require only that there was a segment with good x info.
      // This of course only applies to the Tight algorithms.
      bool existsGoodDTSegX = false;
      bool existsDTSegY = false;

      // Impose cuts on the segments in the station mask until we find a good one
      // Might as well start with the lowest bit to speed things up.
      for(int stationIdx = 0; stationIdx <= 7; ++stationIdx)
         if(theStationMask & 1<<stationIdx) {
            station  = stationIdx < 4 ? stationIdx+1 : stationIdx-3;
            detector = stationIdx < 4 ? 1 : 2;

            if((fabs(muon.pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
                  fabs(muon.dX(station,detector,arbitrationType)) > maxAbsDx) ||
                  (applyAlsoAngularCuts && fabs(muon.pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX))
               continue;
            else if (detector == 1)
               existsGoodDTSegX = true;

            // Is this a tight algorithm?  If yes, use y information
            if (maxAbsDy < 999999) {
               if (detector == 2) { // CSC
                  if((fabs(muon.pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
                        fabs(muon.dY(station,2,arbitrationType)) > maxAbsDy) ||
                        (applyAlsoAngularCuts && fabs(muon.pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY))
                     continue;
               } else {

                  if(muon.dY(station,1,arbitrationType) > 999998) // no y-information
                     continue;
                  else
                     existsDTSegY = true;

                  if((fabs(muon.pullY(station,1,arbitrationType,1)) > maxAbsPullY &&
                        fabs(muon.dY(station,1,arbitrationType)) > maxAbsDy) ||
                        (applyAlsoAngularCuts && fabs(muon.pullDyDz(station,1,arbitrationType,1)) > maxAbsPullY)) {
                     continue;
                  }
               }
            }

            // If we get this far then great this is a good muon
            return true;
         }

      // If we get this far then for sure there are no "good" CSC segments. For
      // DT, check if there were any segments with y information.  If there
      // were none, but there was a segment with good x, then we're happy. If 
      // there WERE segments with y information, then they must have been shit
      // since we are here so fail it.  Of course, if this is a Loose algorithm
      // then fail immediately since if we had good x we would already have
      // returned true
      if (maxAbsDy < 999999) {
         if (existsDTSegY)
            return false;
         else if (existsGoodDTSegX)
            return true;
      } else
         return false;
   } // TMOneStation

   return goodMuon;
}
bool muon::isTightMuon ( const reco::Muon muon,
const reco::Vertex vtx 
)

Definition at line 709 of file MuonSelectors.cc.

References GlobalMuonPromptTight, reco::Muon::innerTrack(), reco::Muon::isGlobalMuon(), isGoodMuon(), reco::Muon::isPFMuon(), reco::Muon::muonBestTrack(), reco::Muon::numberOfMatchedStations(), and reco::Vertex::position().

Referenced by RecoMuonValidator::analyze().

                                                                 {

  if(!muon.isPFMuon() || !muon.isGlobalMuon() ) return false;

  bool muID = isGoodMuon(muon,GlobalMuonPromptTight) && (muon.numberOfMatchedStations() > 1);
    
  
  bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
    muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0; 

  
  bool ip = fabs(muon.muonBestTrack()->dxy(vtx.position())) < 0.2 && fabs(muon.muonBestTrack()->dz(vtx.position())) < 0.5;
  
  return muID && hits && ip;
}
bool muon::overlap ( const reco::Muon muon1,
const reco::Muon muon2,
double  pullX = 1.0,
double  pullY = 1.0,
bool  checkAdjacentChambers = false 
)

Definition at line 640 of file MuonSelectors.cc.

References abs, CSCDetId::chamber(), MuonSubdetId::CSC, CSC(), CSCDetId::endcap(), reco::Muon::matches(), reco::Muon::numberOfMatches(), reco::LeafCandidate::pt(), CSCDetId::ring(), reco::Muon::SegmentAndTrackArbitration, mathSSE::sqrt(), and CSCDetId::station().

Referenced by CSCCFEBData::digis(), TopDiLeptonOffline::MonitorEnsemble::fill(), PFElecTkProducer::isSharingEcalEnergyWithEgSC(), main(), IsolationProducerForTracks::produce(), MatchedProbeMaker< T >::produce(), CSCEventData::selfTest(), CMSMidpointAlgorithm::splitAndMerge(), edm::EventSelector::testSelectionOverlap(), and PFBlockAlgo::testSuperClusterPFCluster().

{
  unsigned int nMatches1 = muon1.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
  unsigned int nMatches2 = muon2.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
  unsigned int betterMuon = ( muon1.pt() > muon2.pt() ? 1 : 2 );
  for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber1 = muon1.matches().begin();
           chamber1 != muon1.matches().end(); ++chamber1 )
    for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber2 = muon2.matches().begin();
          chamber2 != muon2.matches().end(); ++chamber2 )
      {
        
        // if ( (chamber1->segmentMatches.empty() || chamber2->segmentMatches.empty()) ) continue;
        
        // handle case where both muons have information about the same chamber
        // here we know how close they are 
        if ( chamber1->id == chamber2->id ){
          // found the same chamber
          if ( fabs(chamber1->x-chamber2->x) < 
               pullX * sqrt(chamber1->xErr*chamber1->xErr+chamber2->xErr*chamber2->xErr) )
            {
              if ( betterMuon == 1 )
                nMatches2--;
              else
                nMatches1--;
              if ( nMatches1==0 || nMatches2==0 ) return true;
              continue;
            }
          if ( fabs(chamber1->y-chamber2->y) < 
               pullY * sqrt(chamber1->yErr*chamber1->yErr+chamber2->yErr*chamber2->yErr) )
            {
              if ( betterMuon == 1 )
                nMatches2--;
              else
                nMatches1--;
              if ( nMatches1==0 || nMatches2==0 ) return true;
            }
        } else {
          if ( ! checkAdjacentChambers ) continue;
          // check if tracks are pointing into overlaping region of the CSC detector
          if ( chamber1->id.subdetId() != MuonSubdetId::CSC || 
               chamber2->id.subdetId() != MuonSubdetId::CSC ) continue;
          CSCDetId id1(chamber1->id);
          CSCDetId id2(chamber2->id);
          if ( id1.endcap()  != id2.endcap() )  continue;
          if ( id1.station() != id2.station() ) continue;
          if ( id1.ring()    != id2.ring() )    continue;
          if ( abs(id1.chamber() - id2.chamber())>1 ) continue;
          // FIXME: we don't handle 18->1; 36->1 transitions since 
          // I don't know how to check for sure how many chambers
          // are there. Probably need to hard code some checks.
          
          // Now we have to make sure that both tracks are close to an edge
          // FIXME: ignored Y coordinate for now
          if ( fabs(chamber1->edgeX) > chamber1->xErr*pullX ) continue;
          if ( fabs(chamber2->edgeX) > chamber2->xErr*pullX ) continue;
          if ( chamber1->x * chamber2->x < 0 ) { // check if the same edge
            if ( betterMuon == 1 )
              nMatches2--;
            else
              nMatches1--;
            if ( nMatches1==0 || nMatches2==0 ) return true;
          }
        }
      }
  return false;
}
unsigned int muon::RequiredStationMask ( const reco::Muon muon,
double  maxChamberDist,
double  maxChamberDistPull,
reco::Muon::ArbitrationType  arbitrationType 
)

Definition at line 52 of file MuonSelectors.cc.

References reco::Muon::trackDist(), and reco::Muon::trackDistErr().

Referenced by isGoodMuon().

{
   unsigned int theMask = 0;

   for(int stationIdx = 1; stationIdx < 5; ++stationIdx)
      for(int detectorIdx = 1; detectorIdx < 3; ++detectorIdx)
         if(muon.trackDist(stationIdx,detectorIdx,arbitrationType) < maxChamberDist &&
               muon.trackDist(stationIdx,detectorIdx,arbitrationType)/muon.trackDistErr(stationIdx,detectorIdx,arbitrationType) < maxChamberDistPull)
            theMask += 1<<((stationIdx-1)+4*(detectorIdx-1));

   return theMask;
}
float muon::segmentCompatibility ( const reco::Muon muon,
reco::Muon::ArbitrationType  arbitrationType = reco::Muon::SegmentAndTrackArbitration 
)

Definition at line 74 of file MuonSelectors.cc.

References reco::Muon::dX(), reco::Muon::dY(), i, siStripFEDMonitor_P5_cff::Max, reco::Muon::pullX(), reco::Muon::pullY(), reco::Muon::segmentX(), and reco::Muon::trackDist().

Referenced by ExampleMuonAnalyzer::analyze(), MuonIdVal::analyze(), HiggsDQM::analyze(), MuonCosmicCompatibilityFiller::checkMuonSegments(), PFRecoTauDiscriminationAgainstMuon::discriminate(), TauDiscriminationAgainstMuon< TauType, TauDiscriminator >::evaluateMuonVeto(), isGoodMuon(), and pat::MuonSelector::muIdSelection_().

                                                                                             {
  bool use_weight_regain_at_chamber_boundary = true;
  bool use_match_dist_penalty = true;

  int nr_of_stations_crossed = 0;
  int nr_of_stations_with_segment = 0;
  std::vector<int> stations_w_track(8);
  std::vector<int> station_has_segmentmatch(8);
  std::vector<int> station_was_crossed(8);
  std::vector<float> stations_w_track_at_boundary(8);
  std::vector<float> station_weight(8);
  int position_in_stations = 0;
  float full_weight = 0.;

  for(int i = 1; i<=8; ++i) {
    // ********************************************************;
    // *** fill local info for this muon (do some counting) ***;
    // ************** begin ***********************************;
    if(i<=4) { // this is the section for the DTs
      if( muon.trackDist(i,1,arbitrationType) < 999999 ) { //current "raw" info that a track is close to a chamber
        ++nr_of_stations_crossed;
        station_was_crossed[i-1] = 1;
        if(muon.trackDist(i,1,arbitrationType) > -10. ) stations_w_track_at_boundary[i-1] = muon.trackDist(i,1,arbitrationType); 
        else stations_w_track_at_boundary[i-1] = 0.;
      }
      if( muon.segmentX(i,1,arbitrationType) < 999999 ) { //current "raw" info that a segment is matched to the current track
        ++nr_of_stations_with_segment;
        station_has_segmentmatch[i-1] = 1;
      }
    }
    else     { // this is the section for the CSCs
      if( muon.trackDist(i-4,2,arbitrationType) < 999999 ) { //current "raw" info that a track is close to a chamber
        ++nr_of_stations_crossed;
        station_was_crossed[i-1] = 1;
        if(muon.trackDist(i-4,2,arbitrationType) > -10. ) stations_w_track_at_boundary[i-1] = muon.trackDist(i-4,2,arbitrationType);
        else stations_w_track_at_boundary[i-1] = 0.;
      }
      if( muon.segmentX(i-4,2,arbitrationType) < 999999 ) { //current "raw" info that a segment is matched to the current track
        ++nr_of_stations_with_segment;
        station_has_segmentmatch[i-1] = 1;
      }
    }
    // rough estimation of chamber border efficiency (should be parametrized better, this is just a quick guess):
    // TF1 * merf = new TF1("merf","-0.5*(TMath::Erf(x/6.)-1)",-100,100);
    // use above value to "unpunish" missing segment if close to border, i.e. rather than not adding any weight, add
    // the one from the function. Only for dist ~> -10 cm, else full punish!.

    // ********************************************************;
    // *** fill local info for this muon (do some counting) ***;
    // ************** end *************************************;
  }

  // ********************************************************;
  // *** calculate weights for each station *****************;
  // ************** begin ***********************************;
  //    const float slope = 0.5;
  //    const float attenuate_weight_regain = 1.;
  // if attenuate_weight_regain < 1., additional punishment if track is close to boundary and no segment
  const float attenuate_weight_regain = 0.5; 

  for(int i = 1; i<=8; ++i) { // loop over all possible stations

    // first set all weights if a station has been crossed
    // later penalize if a station did not have a matching segment

    //old logic      if(station_has_segmentmatch[i-1] > 0 ) { // the track has an associated segment at the current station
    if( station_was_crossed[i-1] > 0 ) { // the track crossed this chamber (or was nearby)
      // - Apply a weight depending on the "depth" of the muon passage. 
      // - The station_weight is later reduced for stations with badly matched segments. 
      // - Even if there is no segment but the track passes close to a chamber boundary, the
      //   weight is set non zero and can go up to 0.5 of the full weight if the track is quite
      //   far from any station.
      ++position_in_stations;

      switch ( nr_of_stations_crossed ) { // define different weights depending on how many stations were crossed
      case 1 : 
        station_weight[i-1] =  1.;
        break;
      case 2 :
        if     ( position_in_stations == 1 ) station_weight[i-1] =  0.33;
        else                                 station_weight[i-1] =  0.67;
        break;
      case 3 : 
        if     ( position_in_stations == 1 ) station_weight[i-1] =  0.23;
        else if( position_in_stations == 2 ) station_weight[i-1] =  0.33;
        else                                 station_weight[i-1] =  0.44;
        break;
      case 4 : 
        if     ( position_in_stations == 1 ) station_weight[i-1] =  0.10;
        else if( position_in_stations == 2 ) station_weight[i-1] =  0.20;
        else if( position_in_stations == 3 ) station_weight[i-1] =  0.30;
        else                                 station_weight[i-1] =  0.40;
        break;
          
      default : 
//      LogTrace("MuonIdentification")<<"            // Message: A muon candidate track has more than 4 stations with matching segments.";
//      LogTrace("MuonIdentification")<<"            // Did not expect this - please let me know: ibloch@fnal.gov";
        // for all other cases
        station_weight[i-1] = 1./nr_of_stations_crossed;
      }

      if( use_weight_regain_at_chamber_boundary ) { // reconstitute some weight if there is no match but the segment is close to a boundary:
        if(station_has_segmentmatch[i-1] <= 0 && stations_w_track_at_boundary[i-1] != 0. ) {
          // if segment is not present but track in inefficient region, do not count as "missing match" but add some reduced weight. 
          // original "match weight" is currently reduced by at least attenuate_weight_regain, variing with an error function down to 0 if the track is 
          // inside the chamber.
          station_weight[i-1] = station_weight[i-1]*attenuate_weight_regain*0.5*(TMath::Erf(stations_w_track_at_boundary[i-1]/6.)+1.); // remark: the additional scale of 0.5 normalizes Err to run from 0 to 1 in y
        }
        else if(station_has_segmentmatch[i-1] <= 0 && stations_w_track_at_boundary[i-1] == 0.) { // no segment match and track well inside chamber
          // full penalization
          station_weight[i-1] = 0.;
        }
      }
      else { // always fully penalize tracks with no matching segment, whether the segment is close to the boundary or not.
        if(station_has_segmentmatch[i-1] <= 0) station_weight[i-1] = 0.;
      }

      if( station_has_segmentmatch[i-1] > 0 && 42 == 42 ) { // if track has matching segment, but the matching is not high quality, penalize
        if(i<=4) { // we are in the DTs
          if( muon.dY(i,1,arbitrationType) < 999999 && muon.dX(i,1,arbitrationType) < 999999) { // have both X and Y match
            if(
               TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.))> 1. ) {
              // reduce weight
              if(use_match_dist_penalty) {
                // only use pull if 3 sigma is not smaller than 3 cm
                if(TMath::Sqrt(TMath::Power(muon.dX(i,1,arbitrationType),2.)+TMath::Power(muon.dY(i,1,arbitrationType),2.)) < 3. && TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.)) > 3. ) { 
                  station_weight[i-1] *= 1./TMath::Power(
                                                         TMath::Max((double)TMath::Sqrt(TMath::Power(muon.dX(i,1,arbitrationType),2.)+TMath::Power(muon.dY(i,1,arbitrationType),2.)),(double)1.),.25); 
                }
                else {
                  station_weight[i-1] *= 1./TMath::Power(
                                                         TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.)),.25); 
                }
              }
            }
          }
          else if (muon.dY(i,1,arbitrationType) >= 999999) { // has no match in Y
            if( muon.pullX(i,1,arbitrationType) > 1. ) { // has a match in X. Pull larger that 1 to avoid increasing the weight (just penalize, don't anti-penalize)
              // reduce weight
              if(use_match_dist_penalty) {
                // only use pull if 3 sigma is not smaller than 3 cm
                if( muon.dX(i,1,arbitrationType) < 3. && muon.pullX(i,1,arbitrationType) > 3. ) { 
                  station_weight[i-1] *= 1./TMath::Power(TMath::Max((double)muon.dX(i,1,arbitrationType),(double)1.),.25);
                }
                else {
                  station_weight[i-1] *= 1./TMath::Power(muon.pullX(i,1,arbitrationType),.25);
                }
              }
            }
          }
          else { // has no match in X
            if( muon.pullY(i,1,arbitrationType) > 1. ) { // has a match in Y. Pull larger that 1 to avoid increasing the weight (just penalize, don't anti-penalize)
              // reduce weight
              if(use_match_dist_penalty) {
                // only use pull if 3 sigma is not smaller than 3 cm
                if( muon.dY(i,1,arbitrationType) < 3. && muon.pullY(i,1,arbitrationType) > 3. ) { 
                  station_weight[i-1] *= 1./TMath::Power(TMath::Max((double)muon.dY(i,1,arbitrationType),(double)1.),.25);
                }
                else {
                  station_weight[i-1] *= 1./TMath::Power(muon.pullY(i,1,arbitrationType),.25);
                }
              }
            }
          }
        }
        else { // We are in the CSCs
          if(
             TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)) > 1. ) {
            // reduce weight
            if(use_match_dist_penalty) {
              // only use pull if 3 sigma is not smaller than 3 cm
              if(TMath::Sqrt(TMath::Power(muon.dX(i-4,2,arbitrationType),2.)+TMath::Power(muon.dY(i-4,2,arbitrationType),2.)) < 3. && TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)) > 3. ) { 
                station_weight[i-1] *= 1./TMath::Power(
                                                       TMath::Max((double)TMath::Sqrt(TMath::Power(muon.dX(i-4,2,arbitrationType),2.)+TMath::Power(muon.dY(i-4,2,arbitrationType),2.)),(double)1.),.25);
              }
              else {
                station_weight[i-1] *= 1./TMath::Power(
                                                       TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)),.25);
              }
            }
          }
        }
      }
        
      // Thoughts:
      // - should penalize if the segment has only x OR y info
      // - should also use the segment direction, as it now works!
        
    }
    else { // track did not pass a chamber in this station - just reset weight
      station_weight[i-1] = 0.;
    }
      
    //increment final weight for muon:
    full_weight += station_weight[i-1];
  }

  // if we don't expect any matches, we set the compatibility to
  // 0.5 as the track is as compatible with a muon as it is with
  // background - we should maybe rather set it to -0.5!
  if( nr_of_stations_crossed == 0 ) {
    //      full_weight = attenuate_weight_regain*0.5;
    full_weight = 0.5;
  }

  // ********************************************************;
  // *** calculate weights for each station *****************;
  // ************** end *************************************;

  return full_weight;

}
SelectionType muon::selectionTypeFromString ( const std::string &  label)

Definition at line 8 of file MuonSelectors.cc.

References All, AllArbitrated, AllGlobalMuons, AllStandAloneMuons, AllTrackerMuons, Exception, newFWLiteAna::found, GlobalMuonPromptTight, GMStaChiCompatibility, GMTkChiCompatibility, GMTkKinkTight, i, muon::SelectionTypeStringToEnum::label, TM2DCompatibilityLoose, TM2DCompatibilityTight, TMLastStationAngLoose, TMLastStationAngTight, TMLastStationLoose, TMLastStationOptimizedBarrelLowPtLoose, TMLastStationOptimizedBarrelLowPtTight, TMLastStationOptimizedLowPtLoose, TMLastStationOptimizedLowPtTight, TMLastStationTight, TMOneStationAngLoose, TMOneStationAngTight, TMOneStationLoose, TMOneStationTight, TrackerMuonArbitrated, muon::SelectionTypeStringToEnum::value, and relativeConstraints::value.

Referenced by pat::Muon::muonID(), MuonSelectionTypeValueMapProducer::MuonSelectionTypeValueMapProducer(), and MuonTrackProducer::produce().

{
   static SelectionTypeStringToEnum selectionTypeStringToEnumMap[] = {
      { "All", All },
      { "AllGlobalMuons", AllGlobalMuons },
      { "AllStandAloneMuons", AllStandAloneMuons },
      { "AllTrackerMuons", AllTrackerMuons },
      { "TrackerMuonArbitrated", TrackerMuonArbitrated },
      { "AllArbitrated", AllArbitrated },
      { "GlobalMuonPromptTight", GlobalMuonPromptTight },
      { "TMLastStationLoose", TMLastStationLoose },
      { "TMLastStationTight", TMLastStationTight },
      { "TM2DCompatibilityLoose", TM2DCompatibilityLoose },
      { "TM2DCompatibilityTight", TM2DCompatibilityTight },
      { "TMOneStationLoose", TMOneStationLoose },
      { "TMOneStationTight", TMOneStationTight },
      { "TMLastStationOptimizedLowPtLoose", TMLastStationOptimizedLowPtLoose },
      { "TMLastStationOptimizedLowPtTight", TMLastStationOptimizedLowPtTight },
      { "GMTkChiCompatibility", GMTkChiCompatibility },
      { "GMStaChiCompatibility", GMStaChiCompatibility},
      { "GMTkKinkTight", GMTkKinkTight},
      { "TMLastStationAngLoose", TMLastStationAngLoose },
      { "TMLastStationAngTight", TMLastStationAngTight },
      { "TMOneStationAngLoose", TMOneStationAngLoose },
      { "TMOneStationAngTight", TMOneStationAngTight },
      { "TMLastStationOptimizedBarrelLowPtLoose", TMLastStationOptimizedBarrelLowPtLoose },
      { "TMLastStationOptimizedBarrelLowPtTight", TMLastStationOptimizedBarrelLowPtTight },
      { 0, (SelectionType)-1 }
   };

   SelectionType value = (SelectionType)-1;
   bool found = false;
   for(int i = 0; selectionTypeStringToEnumMap[i].label && (! found); ++i)
      if (! strcmp(label.c_str(), selectionTypeStringToEnumMap[i].label)) {
         found = true;
         value = selectionTypeStringToEnumMap[i].value;
      }

   // in case of unrecognized selection type
   if (! found) throw cms::Exception("MuonSelectorError") << label << " is not a recognized SelectionType";
   return value;
}
reco::Muon::MuonTrackTypePair muon::sigmaSwitch ( const reco::TrackRef combinedTrack,
const reco::TrackRef trackerTrack,
const double  nSigma = 2.,
const double  ptThreshold = 200. 
)

Definition at line 93 of file MuonCocktails.cc.

References reco::Muon::CombinedTrack, delta, reco::Muon::InnerTrack, and dtDQMClient_cfg::threshold.

Referenced by sigmaSwitch().

                                                                          {
  // If either the global or tracker-only fits have pT below threshold
  // (default 200 GeV), return the tracker-only fit.
  if (combinedTrack->pt() < ptThreshold || trackerTrack->pt() < ptThreshold)
    return make_pair(trackerTrack,reco::Muon::InnerTrack);
  
  // If both are above the pT threshold, compare the difference in
  // q/p: if less than two sigma of the tracker-only track, switch to
  // global. Otherwise, use tracker-only.
  const double delta = fabs(trackerTrack->qoverp() - combinedTrack->qoverp());
  const double threshold = nSigma * trackerTrack->qoverpError();
  return delta > threshold ? make_pair(trackerTrack,reco::Muon::InnerTrack) :  make_pair(combinedTrack,reco::Muon::CombinedTrack);
}
reco::Muon::MuonTrackTypePair muon::sigmaSwitch ( const reco::Muon muon,
const double  nSigma = 2.,
const double  ptThreshold = 200. 
) [inline]

Definition at line 88 of file MuonCocktails.h.

References reco::Muon::globalTrack(), reco::Muon::innerTrack(), and sigmaSwitch().

                                                                                    {
    return muon::sigmaSwitch(muon.globalTrack(),
                             muon.innerTrack(),
                             nSigma,
                             ptThreshold);
  }
reco::Muon::MuonTrackTypePair muon::tevOptimized ( const reco::Muon muon,
const reco::TrackToTrackMap tevMap1,
const reco::TrackToTrackMap tevMap2,
const reco::TrackToTrackMap tevMap3 
) [inline]

Definition at line 70 of file MuonCocktails.h.

References reco::Muon::combinedMuon(), getTevRefitTrack(), tevOptimized(), and reco::Muon::track().

                                                                                         {
    return tevOptimized(muon.combinedMuon(),
                        muon.track(),
                        getTevRefitTrack(muon.combinedMuon(), tevMap2),
                        getTevRefitTrack(muon.combinedMuon(), tevMap3));
  }
reco::Muon::MuonTrackTypePair muon::tevOptimized ( const reco::TrackRef combinedTrack,
const reco::TrackRef trackerTrack,
const reco::TrackToTrackMap tevMap1,
const reco::TrackToTrackMap tevMap2,
const reco::TrackToTrackMap tevMap3,
const double  ptThreshold = 200.,
const double  tune1 = 4.,
const double  tune2 = 6. 
) [inline]

Definition at line 53 of file MuonCocktails.h.

References getTevRefitTrack(), and tevOptimized().

                                                                             {
    return tevOptimized(combinedTrack,
                        trackerTrack,
                        getTevRefitTrack(combinedTrack, tevMap2),
                        getTevRefitTrack(combinedTrack, tevMap3),
                        ptThreshold,
                        tune1,
                        tune2);
  }
reco::Muon::MuonTrackTypePair muon::tevOptimized ( const reco::Muon muon,
const double  ptThreshold = 200.,
const double  tune1 = 4.,
const double  tune2 = 6. 
) [inline]

Definition at line 30 of file MuonCocktails.h.

References reco::Muon::globalTrack(), reco::Muon::innerTrack(), reco::Muon::pickyTrack(), tevOptimized(), and reco::Muon::tpfmsTrack().

                                                                             {
    return tevOptimized(muon.globalTrack(),
                        muon.innerTrack(),
                        muon.tpfmsTrack(),
                        muon.pickyTrack(),
                        ptThreshold,
                        tune1,
                        tune2);  
  }
reco::Muon::MuonTrackTypePair muon::tevOptimized ( const reco::TrackRef combinedTrack,
const reco::TrackRef trackerTrack,
const reco::TrackRef tpfmsTrack,
const reco::TrackRef pickyTrack,
const double  ptThreshold = 200.,
const double  tune1 = 4.,
const double  tune2 = 6. 
)

Definition at line 9 of file MuonCocktails.cc.

References reco::Muon::CombinedTrack, first, i, reco::Muon::InnerTrack, reco::Muon::Picky, reco::Muon::TPFMS, and trackProbability().

Referenced by MuonIdProducer::makeMuon(), MuonsFromRefitTracksProducer::produce(), and tevOptimized().

                                                                      {
  
  // If Tracker pT is below the pT threshold (currently 200 GeV) - return the Tracker track
  if (trackerTrack->pt() < ptThreshold) return make_pair(trackerTrack,reco::Muon::InnerTrack);  
  
  // Array for convenience below.
  const reco::Muon::MuonTrackTypePair refit[4] = { 
    make_pair(trackerTrack, reco::Muon::InnerTrack), 
    make_pair(combinedTrack,reco::Muon::CombinedTrack),
    make_pair(tpfmsTrack,   reco::Muon::TPFMS),
    make_pair(pickyTrack,   reco::Muon::Picky)
  }; 
  
  // Calculate the log(tail probabilities). If there's a problem,
  // signify this with prob == 0. The current problems recognized are:
  // the track being not available, whether the (re)fit failed or it's
  // just not in the event, or if the (re)fit ended up with no valid
  // hits.
  double prob[4] = {0.};
  for (unsigned int i = 0; i < 4; ++i) 
    if (refit[i].first.isNonnull() && refit[i].first->numberOfValidHits()) 
      prob[i] = muon::trackProbability(refit[i].first); 

  //std::cout << "Probabilities: " << prob[0] << " " << prob[1] << " " << prob[2] << " " << prob[3] << std::endl;
  
  // Start with picky.
  int chosen = 3;
  
  // If there's a problem with picky, make the default one of the
  // other tracks. Try TPFMS first, then global, then tracker-only.
  if (prob[3] == 0.) { 
    if      (prob[2] > 0.) chosen = 2;
    else if (prob[1] > 0.) chosen = 1;
    else if (prob[0] > 0.) chosen = 0;
  } 
  
  // Now the algorithm: switch from picky to tracker-only if the
  // difference, log(tail prob(picky)) - log(tail prob(tracker-only))
  // is greater than a tuned value. Then compare the
  // so-picked track to TPFMS in the same manner using another tuned
  // value.
  if (prob[0] > 0. && prob[3] > 0. && (prob[3] - prob[0]) > tune1)
    chosen = 0;
  if (prob[2] > 0. && (prob[chosen] - prob[2]) > tune2)
    chosen = 2;

  // Done. If pT of the chosen track is below the threshold value, return the tracker track.
  if (refit[chosen].first->pt() < ptThreshold) return make_pair(trackerTrack,reco::Muon::InnerTrack);    
  
  // Return the chosen track (which can be the global track in
  // very rare cases).
  return refit[chosen];
}
reco::Muon::MuonTrackTypePair muon::TMR ( const reco::TrackRef trackerTrack,
const reco::TrackRef fmsTrack,
const double  tune = 4. 
)

Definition at line 113 of file MuonCocktails.cc.

References reco::Muon::InnerTrack, edm::Ref< C, T, F >::isNonnull(), reco::Muon::None, reco::Muon::TPFMS, and trackProbability().

                                                           {
  double probTK  = 0;
  double probFMS = 0;
  
  if (trackerTrack.isNonnull() && trackerTrack->numberOfValidHits())  
    probTK = muon::trackProbability(trackerTrack);
  if (fmsTrack.isNonnull() && fmsTrack->numberOfValidHits())
    probFMS = muon::trackProbability(fmsTrack);
  
  bool TKok  = probTK > 0;
  bool FMSok = probFMS > 0;
  
  if (TKok && FMSok) {
    if (probFMS - probTK > tune)
      return make_pair(trackerTrack,reco::Muon::InnerTrack);
    else
      return  make_pair(fmsTrack,reco::Muon::TPFMS);
  }
  else if (FMSok)
    return  make_pair(fmsTrack,reco::Muon::TPFMS);
  else if (TKok)
    return make_pair(trackerTrack,reco::Muon::InnerTrack);
  else
    return make_pair(reco::TrackRef(),reco::Muon::None); 
}
double muon::trackProbability ( const reco::TrackRef  track)

Definition at line 72 of file MuonCocktails.cc.

References funct::log().

Referenced by HistogramProbabilityEstimator::probability(), tevOptimized(), tevOptimizedTMR(), and TMR().

                                                      {
  
  int nDOF = (int)track->ndof();
  if ( nDOF > 0 && track->chi2()> 0) { 
    return -log(TMath::Prob(track->chi2(), nDOF));
  } else { 
    return 0.0;
  }
  
}