File indexing completed on 2024-04-06 12:26:58
0001
0002
0003 #include "FWCore/Framework/interface/global/EDProducer.h"
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/PluginManager/interface/ModuleDef.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 #include "DataFormats/TrackReco/interface/Track.h"
0012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0013
0014 class L3MuonCleaner : public edm::global::EDProducer<> {
0015 public:
0016 L3MuonCleaner(const edm::ParameterSet&);
0017 ~L3MuonCleaner() override {}
0018 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0019
0020 private:
0021 edm::InputTag m_input;
0022 int m_minTrkHits;
0023 int m_minMuonHits;
0024 double m_maxNormalizedChi2;
0025 edm::EDGetTokenT<reco::TrackCollection> inputToken_;
0026 };
0027
0028 L3MuonCleaner::L3MuonCleaner(const edm::ParameterSet& parameterSet) {
0029 m_input = parameterSet.getParameter<edm::InputTag>("input");
0030 m_minTrkHits = parameterSet.getParameter<int>("minTrkHits");
0031 m_minMuonHits = parameterSet.getParameter<int>("minMuonHits");
0032 m_maxNormalizedChi2 = parameterSet.getParameter<double>("maxNormalizedChi2");
0033 inputToken_ = consumes<reco::TrackCollection>(m_input);
0034
0035 produces<reco::TrackCollection>();
0036 }
0037
0038 void L3MuonCleaner::produce(edm::StreamID, edm::Event& event, const edm::EventSetup&) const {
0039 edm::Handle<reco::TrackCollection> tracks;
0040 event.getByToken(inputToken_, tracks);
0041 auto outTracks = std::make_unique<reco::TrackCollection>();
0042 for (reco::TrackCollection::const_iterator trk = tracks->begin(); trk != tracks->end(); ++trk) {
0043 if (trk->normalizedChi2() > m_maxNormalizedChi2)
0044 continue;
0045 if (trk->hitPattern().numberOfValidTrackerHits() < m_minTrkHits)
0046 continue;
0047 if (trk->hitPattern().numberOfValidMuonHits() < m_minMuonHits)
0048 continue;
0049 outTracks->push_back(*trk);
0050 }
0051 event.put(std::move(outTracks));
0052 }
0053 DEFINE_FWK_MODULE(L3MuonCleaner);