CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/DataFormats/MuonReco/src/MuonSelectors.cc

Go to the documentation of this file.
00001 #include "DataFormats/MuonReco/interface/MuonSelectors.h"
00002 #include "DataFormats/TrackReco/interface/Track.h"
00003 #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
00004 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
00005 #include "DataFormats/VertexReco/interface/Vertex.h"
00006 #include "DataFormats/MuonReco/interface/MuonRPCHitMatch.h"
00007 
00008 namespace muon {
00009 SelectionType selectionTypeFromString( const std::string &label )
00010 {
00011    static SelectionTypeStringToEnum selectionTypeStringToEnumMap[] = {
00012       { "All", All },
00013       { "AllGlobalMuons", AllGlobalMuons },
00014       { "AllStandAloneMuons", AllStandAloneMuons },
00015       { "AllTrackerMuons", AllTrackerMuons },
00016       { "TrackerMuonArbitrated", TrackerMuonArbitrated },
00017       { "AllArbitrated", AllArbitrated },
00018       { "GlobalMuonPromptTight", GlobalMuonPromptTight },
00019       { "TMLastStationLoose", TMLastStationLoose },
00020       { "TMLastStationTight", TMLastStationTight },
00021       { "TM2DCompatibilityLoose", TM2DCompatibilityLoose },
00022       { "TM2DCompatibilityTight", TM2DCompatibilityTight },
00023       { "TMOneStationLoose", TMOneStationLoose },
00024       { "TMOneStationTight", TMOneStationTight },
00025       { "TMLastStationOptimizedLowPtLoose", TMLastStationOptimizedLowPtLoose },
00026       { "TMLastStationOptimizedLowPtTight", TMLastStationOptimizedLowPtTight },
00027       { "GMTkChiCompatibility", GMTkChiCompatibility },
00028       { "GMStaChiCompatibility", GMStaChiCompatibility},
00029       { "GMTkKinkTight", GMTkKinkTight},
00030       { "TMLastStationAngLoose", TMLastStationAngLoose },
00031       { "TMLastStationAngTight", TMLastStationAngTight },
00032       { "TMOneStationAngLoose", TMOneStationAngLoose },
00033       { "TMOneStationAngTight", TMOneStationAngTight },
00034       { "TMLastStationOptimizedBarrelLowPtLoose", TMLastStationOptimizedBarrelLowPtLoose },
00035       { "TMLastStationOptimizedBarrelLowPtTight", TMLastStationOptimizedBarrelLowPtTight },
00036       { "RPCMuLoose", RPCMuLoose },
00037       { 0, (SelectionType)-1 }
00038    };
00039 
00040    SelectionType value = (SelectionType)-1;
00041    bool found = false;
00042    for(int i = 0; selectionTypeStringToEnumMap[i].label && (! found); ++i)
00043       if (! strcmp(label.c_str(), selectionTypeStringToEnumMap[i].label)) {
00044          found = true;
00045          value = selectionTypeStringToEnumMap[i].value;
00046       }
00047 
00048    // in case of unrecognized selection type
00049    if (! found) throw cms::Exception("MuonSelectorError") << label << " is not a recognized SelectionType";
00050    return value;
00051 }
00052 }
00053 
00054 unsigned int muon::RequiredStationMask( const reco::Muon& muon,
00055                                           double maxChamberDist,
00056                                           double maxChamberDistPull,
00057                                           reco::Muon::ArbitrationType arbitrationType )
00058 {
00059    unsigned int theMask = 0;
00060 
00061    for(int stationIdx = 1; stationIdx < 5; ++stationIdx)
00062       for(int detectorIdx = 1; detectorIdx < 3; ++detectorIdx)
00063          if(muon.trackDist(stationIdx,detectorIdx,arbitrationType) < maxChamberDist &&
00064                muon.trackDist(stationIdx,detectorIdx,arbitrationType)/muon.trackDistErr(stationIdx,detectorIdx,arbitrationType) < maxChamberDistPull)
00065             theMask += 1<<((stationIdx-1)+4*(detectorIdx-1));
00066 
00067    return theMask;
00068 }
00069 
00070 // ------------ method to calculate the calo compatibility for a track with matched muon info  ------------
00071 float muon::caloCompatibility(const reco::Muon& muon) {
00072   return muon.caloCompatibility();
00073 }
00074 
00075 // ------------ method to calculate the segment compatibility for a track with matched muon info  ------------
00076 float muon::segmentCompatibility(const reco::Muon& muon,reco::Muon::ArbitrationType arbitrationType) {
00077   bool use_weight_regain_at_chamber_boundary = true;
00078   bool use_match_dist_penalty = true;
00079 
00080   int nr_of_stations_crossed = 0;
00081   int nr_of_stations_with_segment = 0;
00082   std::vector<int> stations_w_track(8);
00083   std::vector<int> station_has_segmentmatch(8);
00084   std::vector<int> station_was_crossed(8);
00085   std::vector<float> stations_w_track_at_boundary(8);
00086   std::vector<float> station_weight(8);
00087   int position_in_stations = 0;
00088   float full_weight = 0.;
00089 
00090   for(int i = 1; i<=8; ++i) {
00091     // ********************************************************;
00092     // *** fill local info for this muon (do some counting) ***;
00093     // ************** begin ***********************************;
00094     if(i<=4) { // this is the section for the DTs
00095       if( muon.trackDist(i,1,arbitrationType) < 999999 ) { //current "raw" info that a track is close to a chamber
00096         ++nr_of_stations_crossed;
00097         station_was_crossed[i-1] = 1;
00098         if(muon.trackDist(i,1,arbitrationType) > -10. ) stations_w_track_at_boundary[i-1] = muon.trackDist(i,1,arbitrationType); 
00099         else stations_w_track_at_boundary[i-1] = 0.;
00100       }
00101       if( muon.segmentX(i,1,arbitrationType) < 999999 ) { //current "raw" info that a segment is matched to the current track
00102         ++nr_of_stations_with_segment;
00103         station_has_segmentmatch[i-1] = 1;
00104       }
00105     }
00106     else     { // this is the section for the CSCs
00107       if( muon.trackDist(i-4,2,arbitrationType) < 999999 ) { //current "raw" info that a track is close to a chamber
00108         ++nr_of_stations_crossed;
00109         station_was_crossed[i-1] = 1;
00110         if(muon.trackDist(i-4,2,arbitrationType) > -10. ) stations_w_track_at_boundary[i-1] = muon.trackDist(i-4,2,arbitrationType);
00111         else stations_w_track_at_boundary[i-1] = 0.;
00112       }
00113       if( muon.segmentX(i-4,2,arbitrationType) < 999999 ) { //current "raw" info that a segment is matched to the current track
00114         ++nr_of_stations_with_segment;
00115         station_has_segmentmatch[i-1] = 1;
00116       }
00117     }
00118     // rough estimation of chamber border efficiency (should be parametrized better, this is just a quick guess):
00119     // TF1 * merf = new TF1("merf","-0.5*(TMath::Erf(x/6.)-1)",-100,100);
00120     // use above value to "unpunish" missing segment if close to border, i.e. rather than not adding any weight, add
00121     // the one from the function. Only for dist ~> -10 cm, else full punish!.
00122 
00123     // ********************************************************;
00124     // *** fill local info for this muon (do some counting) ***;
00125     // ************** end *************************************;
00126   }
00127 
00128   // ********************************************************;
00129   // *** calculate weights for each station *****************;
00130   // ************** begin ***********************************;
00131   //    const float slope = 0.5;
00132   //    const float attenuate_weight_regain = 1.;
00133   // if attenuate_weight_regain < 1., additional punishment if track is close to boundary and no segment
00134   const float attenuate_weight_regain = 0.5; 
00135 
00136   for(int i = 1; i<=8; ++i) { // loop over all possible stations
00137 
00138     // first set all weights if a station has been crossed
00139     // later penalize if a station did not have a matching segment
00140 
00141     //old logic      if(station_has_segmentmatch[i-1] > 0 ) { // the track has an associated segment at the current station
00142     if( station_was_crossed[i-1] > 0 ) { // the track crossed this chamber (or was nearby)
00143       // - Apply a weight depending on the "depth" of the muon passage. 
00144       // - The station_weight is later reduced for stations with badly matched segments. 
00145       // - Even if there is no segment but the track passes close to a chamber boundary, the
00146       //   weight is set non zero and can go up to 0.5 of the full weight if the track is quite
00147       //   far from any station.
00148       ++position_in_stations;
00149 
00150       switch ( nr_of_stations_crossed ) { // define different weights depending on how many stations were crossed
00151       case 1 : 
00152         station_weight[i-1] =  1.;
00153         break;
00154       case 2 :
00155         if     ( position_in_stations == 1 ) station_weight[i-1] =  0.33;
00156         else                                 station_weight[i-1] =  0.67;
00157         break;
00158       case 3 : 
00159         if     ( position_in_stations == 1 ) station_weight[i-1] =  0.23;
00160         else if( position_in_stations == 2 ) station_weight[i-1] =  0.33;
00161         else                                 station_weight[i-1] =  0.44;
00162         break;
00163       case 4 : 
00164         if     ( position_in_stations == 1 ) station_weight[i-1] =  0.10;
00165         else if( position_in_stations == 2 ) station_weight[i-1] =  0.20;
00166         else if( position_in_stations == 3 ) station_weight[i-1] =  0.30;
00167         else                                 station_weight[i-1] =  0.40;
00168         break;
00169           
00170       default : 
00171 //      LogTrace("MuonIdentification")<<"            // Message: A muon candidate track has more than 4 stations with matching segments.";
00172 //      LogTrace("MuonIdentification")<<"            // Did not expect this - please let me know: ibloch@fnal.gov";
00173         // for all other cases
00174         station_weight[i-1] = 1./nr_of_stations_crossed;
00175       }
00176 
00177       if( use_weight_regain_at_chamber_boundary ) { // reconstitute some weight if there is no match but the segment is close to a boundary:
00178         if(station_has_segmentmatch[i-1] <= 0 && stations_w_track_at_boundary[i-1] != 0. ) {
00179           // if segment is not present but track in inefficient region, do not count as "missing match" but add some reduced weight. 
00180           // original "match weight" is currently reduced by at least attenuate_weight_regain, variing with an error function down to 0 if the track is 
00181           // inside the chamber.
00182           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
00183         }
00184         else if(station_has_segmentmatch[i-1] <= 0 && stations_w_track_at_boundary[i-1] == 0.) { // no segment match and track well inside chamber
00185           // full penalization
00186           station_weight[i-1] = 0.;
00187         }
00188       }
00189       else { // always fully penalize tracks with no matching segment, whether the segment is close to the boundary or not.
00190         if(station_has_segmentmatch[i-1] <= 0) station_weight[i-1] = 0.;
00191       }
00192 
00193       if( station_has_segmentmatch[i-1] > 0 && 42 == 42 ) { // if track has matching segment, but the matching is not high quality, penalize
00194         if(i<=4) { // we are in the DTs
00195           if( muon.dY(i,1,arbitrationType) < 999999 && muon.dX(i,1,arbitrationType) < 999999) { // have both X and Y match
00196             if(
00197                TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.))> 1. ) {
00198               // reduce weight
00199               if(use_match_dist_penalty) {
00200                 // only use pull if 3 sigma is not smaller than 3 cm
00201                 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. ) { 
00202                   station_weight[i-1] *= 1./TMath::Power(
00203                                                          TMath::Max((double)TMath::Sqrt(TMath::Power(muon.dX(i,1,arbitrationType),2.)+TMath::Power(muon.dY(i,1,arbitrationType),2.)),(double)1.),.25); 
00204                 }
00205                 else {
00206                   station_weight[i-1] *= 1./TMath::Power(
00207                                                          TMath::Sqrt(TMath::Power(muon.pullX(i,1,arbitrationType),2.)+TMath::Power(muon.pullY(i,1,arbitrationType),2.)),.25); 
00208                 }
00209               }
00210             }
00211           }
00212           else if (muon.dY(i,1,arbitrationType) >= 999999) { // has no match in Y
00213             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)
00214               // reduce weight
00215               if(use_match_dist_penalty) {
00216                 // only use pull if 3 sigma is not smaller than 3 cm
00217                 if( muon.dX(i,1,arbitrationType) < 3. && muon.pullX(i,1,arbitrationType) > 3. ) { 
00218                   station_weight[i-1] *= 1./TMath::Power(TMath::Max((double)muon.dX(i,1,arbitrationType),(double)1.),.25);
00219                 }
00220                 else {
00221                   station_weight[i-1] *= 1./TMath::Power(muon.pullX(i,1,arbitrationType),.25);
00222                 }
00223               }
00224             }
00225           }
00226           else { // has no match in X
00227             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)
00228               // reduce weight
00229               if(use_match_dist_penalty) {
00230                 // only use pull if 3 sigma is not smaller than 3 cm
00231                 if( muon.dY(i,1,arbitrationType) < 3. && muon.pullY(i,1,arbitrationType) > 3. ) { 
00232                   station_weight[i-1] *= 1./TMath::Power(TMath::Max((double)muon.dY(i,1,arbitrationType),(double)1.),.25);
00233                 }
00234                 else {
00235                   station_weight[i-1] *= 1./TMath::Power(muon.pullY(i,1,arbitrationType),.25);
00236                 }
00237               }
00238             }
00239           }
00240         }
00241         else { // We are in the CSCs
00242           if(
00243              TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)) > 1. ) {
00244             // reduce weight
00245             if(use_match_dist_penalty) {
00246               // only use pull if 3 sigma is not smaller than 3 cm
00247               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. ) { 
00248                 station_weight[i-1] *= 1./TMath::Power(
00249                                                        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);
00250               }
00251               else {
00252                 station_weight[i-1] *= 1./TMath::Power(
00253                                                        TMath::Sqrt(TMath::Power(muon.pullX(i-4,2,arbitrationType),2.)+TMath::Power(muon.pullY(i-4,2,arbitrationType),2.)),.25);
00254               }
00255             }
00256           }
00257         }
00258       }
00259         
00260       // Thoughts:
00261       // - should penalize if the segment has only x OR y info
00262       // - should also use the segment direction, as it now works!
00263         
00264     }
00265     else { // track did not pass a chamber in this station - just reset weight
00266       station_weight[i-1] = 0.;
00267     }
00268       
00269     //increment final weight for muon:
00270     full_weight += station_weight[i-1];
00271   }
00272 
00273   // if we don't expect any matches, we set the compatibility to
00274   // 0.5 as the track is as compatible with a muon as it is with
00275   // background - we should maybe rather set it to -0.5!
00276   if( nr_of_stations_crossed == 0 ) {
00277     //      full_weight = attenuate_weight_regain*0.5;
00278     full_weight = 0.5;
00279   }
00280 
00281   // ********************************************************;
00282   // *** calculate weights for each station *****************;
00283   // ************** end *************************************;
00284 
00285   return full_weight;
00286 
00287 }
00288 
00289 bool muon::isGoodMuon( const reco::Muon& muon, 
00290                          AlgorithmType type,
00291                          double minCompatibility,
00292                          reco::Muon::ArbitrationType arbitrationType ) {
00293   if (!muon.isMatchesValid()) return false;
00294   bool goodMuon = false;
00295   
00296   switch( type ) {
00297   case TM2DCompatibility:
00298     // Simplistic first cut in the 2D segment- vs calo-compatibility plane. Will have to be refined!
00299     if( ( (0.8*caloCompatibility( muon ))+(1.2*segmentCompatibility( muon, arbitrationType )) ) > minCompatibility ) goodMuon = true;
00300     else goodMuon = false;
00301     return goodMuon;
00302     break;
00303   default : 
00304     //  LogTrace("MuonIdentification")<<"            // Invalid Algorithm Type called!";
00305     goodMuon = false;
00306     return goodMuon;
00307     break;
00308   }
00309 }
00310 
00311 bool muon::isGoodMuon( const reco::Muon& muon,
00312                          AlgorithmType type,
00313                          int minNumberOfMatches,
00314                          double maxAbsDx,
00315                          double maxAbsPullX,
00316                          double maxAbsDy,
00317                          double maxAbsPullY,
00318                          double maxChamberDist,
00319                          double maxChamberDistPull,
00320                          reco::Muon::ArbitrationType arbitrationType,
00321              bool   syncMinNMatchesNRequiredStationsInBarrelOnly,
00322              bool   applyAlsoAngularCuts)
00323 {
00324    if (!muon.isMatchesValid()) return false;
00325    bool goodMuon = false;
00326 
00327    if (type == TMLastStation) {
00328       // To satisfy my own paranoia, if the user specifies that the
00329       // minimum number of matches is zero, then return true.
00330       if(minNumberOfMatches == 0) return true;
00331 
00332       unsigned int theStationMask = muon.stationMask(arbitrationType);
00333       unsigned int theRequiredStationMask = RequiredStationMask(muon, maxChamberDist, maxChamberDistPull, arbitrationType);
00334 
00335       // Require that there be at least a minimum number of segments
00336       int numSegs = 0;
00337       int numRequiredStations = 0;
00338       for(int it = 0; it < 8; ++it) {
00339          if(theStationMask & 1<<it) ++numSegs;
00340          if(theRequiredStationMask & 1<<it) ++numRequiredStations;
00341       }
00342 
00343       // Make sure the minimum number of matches is not greater than
00344       // the number of required stations but still greater than zero
00345       if (syncMinNMatchesNRequiredStationsInBarrelOnly) {
00346          // Note that we only do this in the barrel region!
00347          if (fabs(muon.eta()) < 1.2) {
00348             if(minNumberOfMatches > numRequiredStations)
00349                minNumberOfMatches = numRequiredStations;
00350             if(minNumberOfMatches < 1) //SK: this only happens for negative values
00351                minNumberOfMatches = 1;
00352          }
00353       } else {
00354          if(minNumberOfMatches > numRequiredStations)
00355             minNumberOfMatches = numRequiredStations;
00356          if(minNumberOfMatches < 1) //SK: this only happens for negative values
00357             minNumberOfMatches = 1;
00358       }
00359 
00360       if(numSegs >= minNumberOfMatches) goodMuon = 1;
00361 
00362       // Require that last required station have segment
00363       // If there are zero required stations keep track
00364       // of the last station with a segment so that we may
00365       // apply the quality cuts below to it instead
00366       int lastSegBit = 0;
00367       if(theRequiredStationMask) {
00368          for(int stationIdx = 7; stationIdx >= 0; --stationIdx)
00369            if(theRequiredStationMask & 1<<stationIdx){
00370                if(theStationMask & 1<<stationIdx) {
00371                   lastSegBit = stationIdx;
00372                   goodMuon &= 1;
00373                   break;
00374                } else {
00375                   goodMuon = false;
00376                   break;
00377                }
00378            }
00379       } else {
00380          for(int stationIdx = 7; stationIdx >= 0; --stationIdx)
00381             if(theStationMask & 1<<stationIdx) {
00382                lastSegBit = stationIdx;
00383                break;
00384             }
00385       }
00386 
00387       if(!goodMuon) return false;
00388 
00389       // Impose pull cuts on last segment
00390       int station = 0, detector = 0;
00391       station  = lastSegBit < 4 ? lastSegBit+1 : lastSegBit-3;
00392       detector = lastSegBit < 4 ? 1 : 2;
00393 
00394       // Check x information
00395       if(fabs(muon.pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
00396             fabs(muon.dX(station,detector,arbitrationType)) > maxAbsDx)
00397          return false;
00398 
00399       if(applyAlsoAngularCuts && fabs(muon.pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX)
00400          return false;
00401 
00402       // Is this a tight algorithm, i.e. do we bother to check y information?
00403       if (maxAbsDy < 999999) { // really if maxAbsDy < 1E9 as currently defined
00404 
00405          // Check y information
00406          if (detector == 2) { // CSC
00407             if(fabs(muon.pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
00408                   fabs(muon.dY(station,2,arbitrationType)) > maxAbsDy)
00409                return false;
00410 
00411             if(applyAlsoAngularCuts && fabs(muon.pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY)
00412                return false;
00413          } else {
00414             //
00415             // In DT, if this is a "Tight" algorithm and the last segment is
00416             // missing y information (always the case in station 4!!!), impose
00417             // respective cuts on the next station in the stationMask that has
00418             // a segment with y information.  If there are no segments with y
00419             // information then there is nothing to penalize. Should we
00420             // penalize in Tight for having zero segments with y information?
00421             // That is the fundamental question.  Of course I am being uber
00422             // paranoid; if this is a good muon then there will probably be at
00423             // least one segment with y information but not always.  Suppose
00424             // somehow a muon only creates segments in station 4, then we
00425             // definitely do not want to require that there be at least one
00426             // segment with y information because we will lose it completely.
00427             //
00428 
00429             for (int stationIdx = station; stationIdx > 0; --stationIdx) {
00430                if(! (theStationMask & 1<<(stationIdx-1)))  // don't bother if the station is not in the stationMask
00431                   continue;
00432 
00433                if(muon.dY(stationIdx,1,arbitrationType) > 999998) // no y-information
00434                   continue;
00435 
00436                if(fabs(muon.pullY(stationIdx,1,arbitrationType,1)) > maxAbsPullY &&
00437                      fabs(muon.dY(stationIdx,1,arbitrationType)) > maxAbsDy) {
00438                   return false;
00439                }
00440 
00441                if(applyAlsoAngularCuts && fabs(muon.pullDyDz(stationIdx,1,arbitrationType,1)) > maxAbsPullY)
00442                   return false;
00443 
00444                // If we get this far then great this is a good muon
00445                return true;
00446             }
00447          }
00448       }
00449 
00450       return goodMuon;
00451    } // TMLastStation
00452 
00453    // TMOneStation requires only that there be one "good" segment, regardless
00454    // of the required stations.  We do not penalize if there are absolutely zero
00455    // segments with y information in the Tight algorithm.  Maybe I'm being
00456    // paranoid but so be it.  If it's really a good muon then we will probably
00457    // find at least one segment with both x and y information but you never
00458    // know, and I don't want to deal with a potential inefficiency in the DT
00459    // like we did with the original TMLastStation.  Incidentally, not penalizing
00460    // for total lack of y information in the Tight algorithm is what is done in
00461    // the new TMLastStation
00462    //                   
00463    if (type == TMOneStation) {
00464       unsigned int theStationMask = muon.stationMask(arbitrationType);
00465 
00466       // Of course there must be at least one segment
00467       if (! theStationMask) return false;
00468 
00469       int  station = 0, detector = 0;
00470       // Keep track of whether or not there is a DT segment with y information.
00471       // In the end, if it turns out there are absolutely zero DT segments with
00472       // y information, require only that there was a segment with good x info.
00473       // This of course only applies to the Tight algorithms.
00474       bool existsGoodDTSegX = false;
00475       bool existsDTSegY = false;
00476 
00477       // Impose cuts on the segments in the station mask until we find a good one
00478       // Might as well start with the lowest bit to speed things up.
00479       for(int stationIdx = 0; stationIdx <= 7; ++stationIdx)
00480          if(theStationMask & 1<<stationIdx) {
00481             station  = stationIdx < 4 ? stationIdx+1 : stationIdx-3;
00482             detector = stationIdx < 4 ? 1 : 2;
00483 
00484             if((fabs(muon.pullX(station,detector,arbitrationType,1)) > maxAbsPullX &&
00485                   fabs(muon.dX(station,detector,arbitrationType)) > maxAbsDx) ||
00486                   (applyAlsoAngularCuts && fabs(muon.pullDxDz(station,detector,arbitrationType,1)) > maxAbsPullX))
00487                continue;
00488             else if (detector == 1)
00489                existsGoodDTSegX = true;
00490 
00491             // Is this a tight algorithm?  If yes, use y information
00492             if (maxAbsDy < 999999) {
00493                if (detector == 2) { // CSC
00494                   if((fabs(muon.pullY(station,2,arbitrationType,1)) > maxAbsPullY &&
00495                         fabs(muon.dY(station,2,arbitrationType)) > maxAbsDy) ||
00496                         (applyAlsoAngularCuts && fabs(muon.pullDyDz(station,2,arbitrationType,1)) > maxAbsPullY))
00497                      continue;
00498                } else {
00499 
00500                   if(muon.dY(station,1,arbitrationType) > 999998) // no y-information
00501                      continue;
00502                   else
00503                      existsDTSegY = true;
00504 
00505                   if((fabs(muon.pullY(station,1,arbitrationType,1)) > maxAbsPullY &&
00506                         fabs(muon.dY(station,1,arbitrationType)) > maxAbsDy) ||
00507                         (applyAlsoAngularCuts && fabs(muon.pullDyDz(station,1,arbitrationType,1)) > maxAbsPullY)) {
00508                      continue;
00509                   }
00510                }
00511             }
00512 
00513             // If we get this far then great this is a good muon
00514             return true;
00515          }
00516 
00517       // If we get this far then for sure there are no "good" CSC segments. For
00518       // DT, check if there were any segments with y information.  If there
00519       // were none, but there was a segment with good x, then we're happy. If 
00520       // there WERE segments with y information, then they must have been shit
00521       // since we are here so fail it.  Of course, if this is a Loose algorithm
00522       // then fail immediately since if we had good x we would already have
00523       // returned true
00524       if (maxAbsDy < 999999) {
00525          if (existsDTSegY)
00526             return false;
00527          else if (existsGoodDTSegX)
00528             return true;
00529       } else
00530          return false;
00531    } // TMOneStation
00532 
00533   if ( type == RPCMu )
00534   {
00535     if ( minNumberOfMatches == 0 ) return true;
00536     
00537     int nMatch = 0;
00538     for ( std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = muon.matches().begin();
00539           chamberMatch != muon.matches().end(); ++chamberMatch )
00540     {
00541       if ( chamberMatch->detector() != 3 ) continue;
00542 
00543       const double trkX = chamberMatch->x;
00544       const double errX = chamberMatch->xErr;
00545 
00546       for ( std::vector<reco::MuonRPCHitMatch>::const_iterator rpcMatch = chamberMatch->rpcMatches.begin();
00547             rpcMatch != chamberMatch->rpcMatches.end(); ++rpcMatch )
00548       {
00549         const double rpcX = rpcMatch->x;
00550 
00551         const double dX = fabs(rpcX-trkX);
00552         if ( dX < maxAbsDx or dX/errX < maxAbsPullX )
00553         {
00554           ++nMatch;
00555           break;
00556         }
00557       }
00558     }
00559 
00560     if ( nMatch >= minNumberOfMatches ) return true;
00561     else return false;
00562   } // RPCMu
00563 
00564    return goodMuon;
00565 }
00566 
00567 bool muon::isGoodMuon( const reco::Muon& muon, SelectionType type,
00568                        reco::Muon::ArbitrationType arbitrationType)
00569 {
00570   switch (type)
00571     {
00572     case muon::All:
00573       return true;
00574       break;
00575     case muon::AllGlobalMuons:
00576       return muon.isGlobalMuon();
00577       break;
00578     case muon::AllTrackerMuons:
00579       return muon.isTrackerMuon();
00580       break;
00581     case muon::AllStandAloneMuons:
00582       return muon.isStandAloneMuon();
00583       break;
00584     case muon::TrackerMuonArbitrated:
00585       return muon.isTrackerMuon() && muon.numberOfMatches(arbitrationType)>0;
00586       break;
00587     case muon::AllArbitrated:
00588       return ! muon.isTrackerMuon() || muon.numberOfMatches(arbitrationType)>0;
00589       break;
00590     case muon::GlobalMuonPromptTight:
00591       return muon.isGlobalMuon() && muon.globalTrack()->normalizedChi2()<10. && muon.globalTrack()->hitPattern().numberOfValidMuonHits() >0;
00592       break;
00593       // For "Loose" algorithms we choose maximum y quantity cuts of 1E9 instead of
00594       // 9999 as before.  We do this because the muon methods return 999999 (note
00595       // there are six 9's) when the requested information is not available.  For
00596       // example, if a muon fails to traverse the z measuring superlayer in a station
00597       // in the DT, then all methods involving segmentY in this station return
00598       // 999999 to demonstrate that the information is missing.  In order to not
00599       // penalize muons for missing y information in Loose algorithms where we do
00600       // not care at all about y information, we raise these limits.  In the
00601       // TMLastStation and TMOneStation algorithms we actually use this huge number
00602       // to determine whether to consider y information at all.
00603     case muon::TMLastStationLoose:
00604       return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,true,false);
00605       break;
00606     case muon::TMLastStationTight:
00607       return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,true,false);
00608       break;
00609     case muon::TMOneStationLoose:
00610       return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
00611       break;
00612     case muon::TMOneStationTight:
00613       return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
00614       break;
00615     case muon::TMLastStationOptimizedLowPtLoose:
00616       if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
00617         return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
00618       else
00619         return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,false,false);
00620       break;
00621     case muon::TMLastStationOptimizedLowPtTight:
00622       if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
00623         return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
00624       else
00625         return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,false,false);
00626       break;
00627       //compatibility loose
00628     case muon::TM2DCompatibilityLoose:
00629       return muon.isTrackerMuon() && isGoodMuon(muon,TM2DCompatibility,0.7,arbitrationType);
00630       break;
00631       //compatibility tight
00632     case muon::TM2DCompatibilityTight:
00633       return muon.isTrackerMuon() && isGoodMuon(muon,TM2DCompatibility,1.0,arbitrationType);
00634       break;
00635     case muon::GMTkChiCompatibility:
00636       return muon.isGlobalMuon() && muon.isQualityValid() && fabs(muon.combinedQuality().trkRelChi2 - muon.innerTrack()->normalizedChi2()) < 2.0;
00637       break;
00638     case muon::GMStaChiCompatibility:
00639       return muon.isGlobalMuon() && muon.isQualityValid() && fabs(muon.combinedQuality().staRelChi2 - muon.outerTrack()->normalizedChi2()) < 2.0;
00640       break;
00641     case muon::GMTkKinkTight:
00642       return muon.isGlobalMuon() && muon.isQualityValid() && muon.combinedQuality().trkKink < 100.0;
00643       break;
00644     case muon::TMLastStationAngLoose:
00645       return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,false,true);
00646       break;
00647     case muon::TMLastStationAngTight:
00648       return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,false,true);
00649       break;
00650     case muon::TMOneStationAngLoose:
00651       return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,true);
00652       break;
00653     case muon::TMOneStationAngTight:
00654       return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,true);
00655       break;
00656     case muon::TMLastStationOptimizedBarrelLowPtLoose:
00657       if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
00658         return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,1E9,1E9,1E9,1E9,arbitrationType,false,false);
00659       else
00660         return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,1E9,1E9,-3,-3,arbitrationType,true,false);
00661       break;
00662     case muon::TMLastStationOptimizedBarrelLowPtTight:
00663       if (muon.pt() < 8. && fabs(muon.eta()) < 1.2)
00664         return muon.isTrackerMuon() && isGoodMuon(muon,TMOneStation,1,3,3,3,3,1E9,1E9,arbitrationType,false,false);
00665       else
00666         return muon.isTrackerMuon() && isGoodMuon(muon,TMLastStation,2,3,3,3,3,-3,-3,arbitrationType,true,false);
00667       break;
00668     case muon::RPCMuLoose:
00669         return muon.isRPCMuon() && isGoodMuon(muon, RPCMu, 2, 20, 4, 1e9, 1e9, 1e9, 1e9, arbitrationType, false, false);
00670       break;
00671     default:
00672       return false;
00673     }
00674 }
00675 
00676 bool muon::overlap( const reco::Muon& muon1, const reco::Muon& muon2, 
00677                     double pullX, double pullY, bool checkAdjacentChambers)
00678 {
00679   unsigned int nMatches1 = muon1.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
00680   unsigned int nMatches2 = muon2.numberOfMatches(reco::Muon::SegmentAndTrackArbitration);
00681   unsigned int betterMuon = ( muon1.pt() > muon2.pt() ? 1 : 2 );
00682   for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber1 = muon1.matches().begin();
00683            chamber1 != muon1.matches().end(); ++chamber1 )
00684     for ( std::vector<reco::MuonChamberMatch>::const_iterator chamber2 = muon2.matches().begin();
00685           chamber2 != muon2.matches().end(); ++chamber2 )
00686       {
00687         
00688         // if ( (chamber1->segmentMatches.empty() || chamber2->segmentMatches.empty()) ) continue;
00689         
00690         // handle case where both muons have information about the same chamber
00691         // here we know how close they are 
00692         if ( chamber1->id == chamber2->id ){
00693           // found the same chamber
00694           if ( fabs(chamber1->x-chamber2->x) < 
00695                pullX * sqrt(chamber1->xErr*chamber1->xErr+chamber2->xErr*chamber2->xErr) )
00696             {
00697               if ( betterMuon == 1 )
00698                 nMatches2--;
00699               else
00700                 nMatches1--;
00701               if ( nMatches1==0 || nMatches2==0 ) return true;
00702               continue;
00703             }
00704           if ( fabs(chamber1->y-chamber2->y) < 
00705                pullY * sqrt(chamber1->yErr*chamber1->yErr+chamber2->yErr*chamber2->yErr) )
00706             {
00707               if ( betterMuon == 1 )
00708                 nMatches2--;
00709               else
00710                 nMatches1--;
00711               if ( nMatches1==0 || nMatches2==0 ) return true;
00712             }
00713         } else {
00714           if ( ! checkAdjacentChambers ) continue;
00715           // check if tracks are pointing into overlaping region of the CSC detector
00716           if ( chamber1->id.subdetId() != MuonSubdetId::CSC || 
00717                chamber2->id.subdetId() != MuonSubdetId::CSC ) continue;
00718           CSCDetId id1(chamber1->id);
00719           CSCDetId id2(chamber2->id);
00720           if ( id1.endcap()  != id2.endcap() )  continue;
00721           if ( id1.station() != id2.station() ) continue;
00722           if ( id1.ring()    != id2.ring() )    continue;
00723           if ( abs(id1.chamber() - id2.chamber())>1 ) continue;
00724           // FIXME: we don't handle 18->1; 36->1 transitions since 
00725           // I don't know how to check for sure how many chambers
00726           // are there. Probably need to hard code some checks.
00727           
00728           // Now we have to make sure that both tracks are close to an edge
00729           // FIXME: ignored Y coordinate for now
00730           if ( fabs(chamber1->edgeX) > chamber1->xErr*pullX ) continue;
00731           if ( fabs(chamber2->edgeX) > chamber2->xErr*pullX ) continue;
00732           if ( chamber1->x * chamber2->x < 0 ) { // check if the same edge
00733             if ( betterMuon == 1 )
00734               nMatches2--;
00735             else
00736               nMatches1--;
00737             if ( nMatches1==0 || nMatches2==0 ) return true;
00738           }
00739         }
00740       }
00741   return false;
00742 }
00743 
00744 
00745 bool muon::isTightMuon(const reco::Muon& muon, const reco::Vertex& vtx){
00746 
00747   if(!muon.isPFMuon() || !muon.isGlobalMuon() ) return false;
00748 
00749   bool muID = isGoodMuon(muon,GlobalMuonPromptTight) && (muon.numberOfMatchedStations() > 1);
00750     
00751   
00752   bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
00753     muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0; 
00754 
00755   
00756   bool ip = fabs(muon.muonBestTrack()->dxy(vtx.position())) < 0.2 && fabs(muon.muonBestTrack()->dz(vtx.position())) < 0.5;
00757   
00758   return muID && hits && ip;
00759 }
00760 
00761 
00762 bool muon::isLooseMuon(const reco::Muon& muon){
00763   return muon.isPFMuon() && ( muon.isGlobalMuon() || muon.isTrackerMuon());
00764 }
00765 
00766 
00767 bool muon::isSoftMuon(const reco::Muon& muon, const reco::Vertex& vtx){
00768 
00769   bool muID = muon::isGoodMuon(muon, TMOneStationTight);
00770 
00771   if(!muID) return false;
00772   
00773   bool layers = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
00774     muon.innerTrack()->hitPattern().pixelLayersWithMeasurement() > 1;
00775 
00776   bool chi2 = muon.innerTrack()->normalizedChi2() < 1.8;  
00777   
00778   bool ip = fabs(muon.innerTrack()->dxy(vtx.position())) < 3. && fabs(muon.innerTrack()->dz(vtx.position())) < 30.;
00779   
00780   return muID && layers && ip && chi2 ;
00781 }
00782 
00783 
00784 
00785 bool muon::isHighPtMuon(const reco::Muon& muon, const reco::Vertex& vtx){
00786   bool muID =   muon.isGlobalMuon() && muon.globalTrack()->hitPattern().numberOfValidMuonHits() >0 && (muon.numberOfMatchedStations() > 1);
00787   if(!muID) return false;
00788 
00789   bool hits = muon.innerTrack()->hitPattern().trackerLayersWithMeasurement() > 5 &&
00790     muon.innerTrack()->hitPattern().numberOfValidPixelHits() > 0; 
00791 
00792   bool momQuality = muon.muonBestTrack()->ptError()/muon.muonBestTrack()->pt() < 0.3;
00793 
00794   bool ip = fabs(muon.muonBestTrack()->dxy(vtx.position())) < 0.2 && fabs(muon.bestTrack()->dz(vtx.position())) < 0.5;
00795   
00796   return muID && hits && momQuality && ip;
00797 
00798 }
00799 
00800 int muon::sharedSegments( const reco::Muon& mu, const reco::Muon& mu2, unsigned int segmentArbitrationMask ) {
00801     int ret = 0;
00802    
00803     // Will do with a stupid double loop, since creating and filling a map is probably _more_ inefficient for a single lookup.
00804     for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch = mu.matches().begin();
00805         chamberMatch != mu.matches().end(); ++chamberMatch) {
00806         if (chamberMatch->segmentMatches.empty()) continue;
00807         for(std::vector<reco::MuonChamberMatch>::const_iterator chamberMatch2 = mu2.matches().begin();
00808             chamberMatch2 != mu2.matches().end(); ++chamberMatch2) {
00809             if (chamberMatch2->segmentMatches.empty()) continue;
00810             if (chamberMatch2->id() != chamberMatch->id()) continue;
00811             for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch = chamberMatch->segmentMatches.begin(); 
00812                 segmentMatch != chamberMatch->segmentMatches.end(); ++segmentMatch) {
00813                 if (!segmentMatch->isMask(segmentArbitrationMask)) continue;
00814                 for(std::vector<reco::MuonSegmentMatch>::const_iterator segmentMatch2 = chamberMatch2->segmentMatches.begin(); 
00815                     segmentMatch2 != chamberMatch2->segmentMatches.end(); ++segmentMatch2) {
00816                     if (!segmentMatch2->isMask(segmentArbitrationMask)) continue;
00817                     if ((segmentMatch->cscSegmentRef.isNonnull() && segmentMatch->cscSegmentRef == segmentMatch2->cscSegmentRef) ||
00818                         (segmentMatch-> dtSegmentRef.isNonnull() && segmentMatch-> dtSegmentRef == segmentMatch2-> dtSegmentRef) ) {
00819                         ++ret;
00820                     } // is the same
00821                 } // segment of mu2 in chamber
00822             } // segment of mu1 in chamber
00823         } // chamber of mu2
00824     } // chamber of mu1
00825   
00826     return ret; 
00827 }