CMS 3D CMS Logo

HLTMuonTrackSelector.cc
Go to the documentation of this file.
1 /*
2 * class HLTMuonTrackSelector
3 *
4 * See header file for documentation
5 *
6 * Author: Kyeongpil Lee (kplee@cern.ch)
7 *
8 */
9 #include <cmath>
10 
11 #include "HLTMuonTrackSelector.h"
12 
18 
19 using namespace std;
20 using namespace reco;
21 
23  : collectionCloner(producesCollector(), iConfig, true),
24  collectionClonerTokens(iConfig.getParameter<edm::InputTag>("track"), consumesCollector()),
25  token_muon(consumes<vector<reco::Muon> >(iConfig.getParameter<edm::InputTag>("muon"))),
26  token_originalMVAVals(consumes<MVACollection>(iConfig.getParameter<edm::InputTag>("originalMVAVals"))),
27  flag_copyMVA(iConfig.getParameter<bool>("copyMVA")) {
28  produces<MVACollection>("MVAValues");
29 }
30 
33  desc.add<edm::InputTag>("track", edm::InputTag());
34  desc.add<edm::InputTag>("muon", edm::InputTag());
35  desc.add<edm::InputTag>("originalMVAVals", edm::InputTag());
36  desc.add<bool>("copyMVA", false);
37  TrackCollectionCloner::fill(desc); // -- add copyExtras and copyTrajectories
38  descriptions.add("HLTMuonTrackSelector", desc);
39 }
40 
43 
44  // -- load tracks
46  auto nTrack = originalTracks.size();
47 
48  // -- load muons
49  edm::Handle<vector<reco::Muon> > handle_muon;
50  iEvent.getByToken(token_muon, handle_muon);
51  auto nMuon = handle_muon->size();
52 
53  // -- load MVA values if necessary
54  edm::Handle<MVACollection> handle_originalMVAVals;
55  if (flag_copyMVA) {
56  iEvent.getByToken(token_originalMVAVals, handle_originalMVAVals);
57  assert((*handle_originalMVAVals).size() == nTrack);
58  }
59 
60  // -- containers for selected track informations
61  std::vector<unsigned int> selectedIter;
62  auto selectedMVAVals = std::make_unique<MVACollection>();
63 
64  auto nSelected = 0U;
65 
67  // -- matching -- //
69 
70  // -- iteration over muons
71  for (auto i_mu = 0U; i_mu < nMuon; ++i_mu) {
72  // -- avoids crashing in case the muon is SA only.
73  const reco::Muon& muon(handle_muon->at(i_mu));
74  TrackRef muonTrackRef = (muon.innerTrack().isNonnull()) ? muon.innerTrack() : muon.muonBestTrack();
75 
76  double muonPt = muonTrackRef->pt();
77  double muonEta = muonTrackRef->eta();
78  double muonPhi = muonTrackRef->phi();
79 
80  double smallestDPt = 1e30;
81  unsigned int smallestDPtIter = 9999U;
82 
83  // -- iteration over tracks
84  for (auto i_trk = 0U; i_trk < nTrack; ++i_trk) {
85  auto const& track = originalTracks[i_trk];
86 
87  double trackPt = track.pt();
88  double trackEta = track.eta();
89  double trackPhi = track.phi();
90 
92  double dPt = fabs(trackPt - muonPt);
93 
94  if (dPt < smallestDPt) {
95  smallestDPt = dPt;
96  smallestDPtIter = i_trk;
97  }
98  }
99  } // -- end of track iteration
100 
101  // -- if at least one track is matched
102  if (smallestDPtIter != 9999U) {
103  selectedIter.push_back(smallestDPtIter);
104  if (flag_copyMVA)
105  selectedMVAVals->push_back((*handle_originalMVAVals)[smallestDPtIter]);
106  ++nSelected;
107  }
108 
109  } // -- end of muon iteration
110 
111  assert(producer.selTracks_->empty());
112 
113  // -- produces tracks and associated informations
114  producer(collectionClonerTokens, selectedIter);
115  assert(producer.selTracks_->size() == nSelected);
116 
117  if (flag_copyMVA)
118  iEvent.put(std::move(selectedMVAVals), "MVAValues");
119 }
120 
122 
const TrackCollectionCloner::Tokens collectionClonerTokens
const edm::EDGetTokenT< std::vector< reco::Muon > > token_muon
static void fill(edm::ParameterSetDescription &desc)
assert(be >=bs)
static void fillDescriptions(edm::ConfigurationDescriptions &)
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< MVACollection > token_originalMVAVals
Definition: Muon.py:1
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HLTMuonTrackSelector(const edm::ParameterSet &)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
originalTracks
Definition: MTS_cfg.py:184
reco::TrackCollection const & tracks(edm::Event &evt) const
def move(src, dest)
Definition: eostools.py:511
TrackCollectionCloner collectionCloner
std::vector< float > MVACollection