Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:44

0001 #ifndef __RecoParticleFlow_Benchmark_Matchers__
0002 #define __RecoParticleFlow_Benchmark_Matchers__
0003 
0004 #include "DataFormats/Math/interface/deltaR.h"
0005 
0006 #include "DataFormats/MuonReco/interface/MuonSelectors.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include "FWCore/Utilities/interface/EDGetToken.h"
0010 
0011 #include <iostream>
0012 #include <vector>
0013 
0014 namespace PFB {
0015 
0016   template <typename C, typename M>
0017   void match(const C &candCollection,
0018              const M &matchedCandCollection,
0019              std::vector<int> &matchIndices,
0020              bool matchCharge = false,
0021              float dRMax = -1) {
0022     // compute distance to each candidate in the matchedCandCollection.
0023 
0024     float dR2Max = 0;
0025     if (dRMax > 0)
0026       dR2Max = dRMax * dRMax;
0027 
0028     matchIndices.clear();
0029     matchIndices.resize(candCollection.size(), -1);
0030 
0031     for (unsigned i = 0; i < candCollection.size(); ++i) {
0032       static const double bigNumber = 1e14;
0033       double dR2min = bigNumber;
0034       int jMin = -1;
0035       for (unsigned jm = 0; jm < matchedCandCollection.size(); ++jm) {
0036         if (matchCharge && candCollection[i].charge() != matchedCandCollection[jm].charge())
0037           continue;
0038 
0039         double dR2 = reco::deltaR2(candCollection[i], matchedCandCollection[jm]);
0040         if (dR2 < dR2min) {
0041           dR2min = dR2;
0042           jMin = jm;
0043         }
0044       }
0045 
0046       if ((dR2Max > 0 && dR2min < dR2Max) || dRMax <= 0) {
0047         matchIndices[i] = jMin;
0048         /*  std::cout<<"match "<<dR2min<<std::endl;  */
0049       }
0050       // store the closest match, no cut on deltaR.
0051     }
0052   }
0053 
0054   // needed by muon matching
0055   // template< typename C, typename M, typename MM>
0056   template <typename C, typename M>
0057   void match(const C &candCollection,
0058              const M &matchedCandCollection,
0059              std::vector<int> &matchIndices,
0060              const edm::ParameterSet &parameterSet,
0061              // const MM& muonMatchedCandCollection,
0062              edm::View<reco::Muon> muonMatchedCandCollection,
0063              bool matchCharge = false,
0064              float dRMax = -1) {
0065     // compute distance to each candidate in the matchedCandCollection.
0066 
0067     float dR2Max = 0;
0068     if (dRMax > 0)
0069       dR2Max = dRMax * dRMax;
0070 
0071     matchIndices.clear();
0072     matchIndices.resize(candCollection.size(), -1);
0073 
0074     for (unsigned i = 0; i < candCollection.size(); ++i) {
0075       static const double bigNumber = 1e14;
0076       double dR2min = bigNumber;
0077       int jMin = -1;
0078       for (unsigned jm = 0; jm < matchedCandCollection.size(); ++jm) {
0079         if (parameterSet.getParameter<bool>("slimmedLikeSelection")) {
0080           if (!(muonMatchedCandCollection[jm].pt() > parameterSet.getParameter<double>("ptBase") ||
0081                 muonMatchedCandCollection[jm].isPFMuon() ||
0082                 (muonMatchedCandCollection[jm].pt() > parameterSet.getParameter<double>("ptNotPF") &&
0083                  (muonMatchedCandCollection[jm].isGlobalMuon() || muonMatchedCandCollection[jm].isStandAloneMuon() ||
0084                   muonMatchedCandCollection[jm].numberOfMatches() > 0 ||
0085                   muon::isGoodMuon(muonMatchedCandCollection[jm], muon::RPCMuLoose)))))
0086             continue;
0087         }
0088 
0089         if (matchCharge && candCollection[i].charge() != matchedCandCollection[jm].charge())
0090           continue;
0091 
0092         double dR2 = reco::deltaR2(candCollection[i], matchedCandCollection[jm]);
0093         if (dR2 < dR2min) {
0094           dR2min = dR2;
0095           jMin = jm;
0096         }
0097       }
0098 
0099       if ((dR2Max > 0 && dR2min < dR2Max) || dRMax <= 0) {
0100         matchIndices[i] = jMin;
0101         /*  std::cout<<"match "<<dR2min<<std::endl;  */
0102       }
0103       // store the closest match, no cut on deltaR.
0104     }
0105   }
0106 
0107 }  // namespace PFB
0108 
0109 #endif