Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-09 23:33:53

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