Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:38

0001 /*
0002 * class HLTMuonTrackSelector
0003 * 
0004 * See header file for documentation
0005 *  
0006 * Author: Kyeongpil Lee (kplee@cern.ch)
0007 * 
0008 */
0009 #include <cmath>
0010 
0011 #include "HLTMuonTrackSelector.h"
0012 
0013 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0014 #include "DataFormats/TrackReco/interface/Track.h"
0015 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0016 #include "DataFormats/MuonReco/interface/Muon.h"
0017 #include "DataFormats/Math/interface/deltaR.h"
0018 
0019 using namespace std;
0020 using namespace reco;
0021 
0022 HLTMuonTrackSelector::HLTMuonTrackSelector(const edm::ParameterSet& iConfig)
0023     : collectionCloner(producesCollector(), iConfig, true),
0024       collectionClonerTokens(iConfig.getParameter<edm::InputTag>("track"), consumesCollector()),
0025       token_muon(consumes<vector<reco::Muon> >(iConfig.getParameter<edm::InputTag>("muon"))),
0026       token_originalMVAVals(consumes<MVACollection>(iConfig.getParameter<edm::InputTag>("originalMVAVals"))),
0027       flag_copyMVA(iConfig.getParameter<bool>("copyMVA")) {
0028   produces<MVACollection>("MVAValues");
0029 }
0030 
0031 void HLTMuonTrackSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0032   edm::ParameterSetDescription desc;
0033   desc.add<edm::InputTag>("track", edm::InputTag());
0034   desc.add<edm::InputTag>("muon", edm::InputTag());
0035   desc.add<edm::InputTag>("originalMVAVals", edm::InputTag());
0036   desc.add<bool>("copyMVA", false);
0037   TrackCollectionCloner::fill(desc);  // -- add copyExtras and copyTrajectories
0038   descriptions.add("HLTMuonTrackSelector", desc);
0039 }
0040 
0041 void HLTMuonTrackSelector::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0042   TrackCollectionCloner::Producer producer(iEvent, collectionCloner);
0043 
0044   // -- load tracks
0045   auto const& originalTracks = collectionClonerTokens.tracks(iEvent);
0046   auto nTrack = originalTracks.size();
0047 
0048   // -- load muons
0049   edm::Handle<vector<reco::Muon> > handle_muon;
0050   iEvent.getByToken(token_muon, handle_muon);
0051   auto nMuon = handle_muon->size();
0052 
0053   // -- load MVA values if necessary
0054   edm::Handle<MVACollection> handle_originalMVAVals;
0055   if (flag_copyMVA) {
0056     iEvent.getByToken(token_originalMVAVals, handle_originalMVAVals);
0057     assert((*handle_originalMVAVals).size() == nTrack);
0058   }
0059 
0060   // -- containers for selected track informations
0061   std::vector<unsigned int> selectedIter;
0062   auto selectedMVAVals = std::make_unique<MVACollection>();
0063 
0064   auto nSelected = 0U;
0065 
0066   ////////////////////
0067   // -- matching -- //
0068   ////////////////////
0069 
0070   // -- iteration over muons
0071   for (auto i_mu = 0U; i_mu < nMuon; ++i_mu) {
0072     // -- avoids crashing in case the muon is SA only.
0073     const reco::Muon& muon(handle_muon->at(i_mu));
0074     TrackRef muonTrackRef = (muon.innerTrack().isNonnull()) ? muon.innerTrack() : muon.muonBestTrack();
0075 
0076     double muonPt = muonTrackRef->pt();
0077     double muonEta = muonTrackRef->eta();
0078     double muonPhi = muonTrackRef->phi();
0079 
0080     double smallestDPt = 1e30;
0081     unsigned int smallestDPtIter = 9999U;
0082 
0083     // -- iteration over tracks
0084     for (auto i_trk = 0U; i_trk < nTrack; ++i_trk) {
0085       auto const& track = originalTracks[i_trk];
0086 
0087       double trackPt = track.pt();
0088       double trackEta = track.eta();
0089       double trackPhi = track.phi();
0090 
0091       if (reco::deltaR2(trackEta, trackPhi, muonEta, muonPhi) < 0.01) {
0092         double dPt = fabs(trackPt - muonPt);
0093 
0094         if (dPt < smallestDPt) {
0095           smallestDPt = dPt;
0096           smallestDPtIter = i_trk;
0097         }
0098       }
0099     }  // -- end of track iteration
0100 
0101     // -- if at least one track is matched
0102     if (smallestDPtIter != 9999U) {
0103       selectedIter.push_back(smallestDPtIter);
0104       if (flag_copyMVA)
0105         selectedMVAVals->push_back((*handle_originalMVAVals)[smallestDPtIter]);
0106       ++nSelected;
0107     }
0108 
0109   }  // -- end of muon iteration
0110 
0111   assert(producer.selTracks_->empty());
0112 
0113   // -- produces tracks and associated informations
0114   producer(collectionClonerTokens, selectedIter);
0115   assert(producer.selTracks_->size() == nSelected);
0116 
0117   if (flag_copyMVA)
0118     iEvent.put(std::move(selectedMVAVals), "MVAValues");
0119 }
0120 
0121 #include "FWCore/Framework/interface/MakerMacros.h"
0122 
0123 DEFINE_FWK_MODULE(HLTMuonTrackSelector);