CMS 3D CMS Logo

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