File indexing completed on 2024-04-06 12:25:28
0001
0002
0003
0004
0005 #include "FWCore/Framework/interface/global/EDFilter.h"
0006 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 #include <atomic>
0010
0011 class ptHatFilter : public edm::global::EDFilter<> {
0012 public:
0013 ptHatFilter(const edm::ParameterSet&);
0014 void beginJob() override;
0015 bool filter(edm::StreamID, edm::Event& e, edm::EventSetup const& iSetup) const override;
0016 void endJob() override;
0017
0018 private:
0019 const edm::EDGetTokenT<double> token_;
0020 double ptHatLowerCut;
0021 double ptHatUpperCut;
0022 mutable std::atomic<int> totalEvents;
0023 mutable std::atomic<int> acceptedEvents;
0024 };
0025
0026 using namespace edm;
0027 using namespace reco;
0028 using namespace std;
0029
0030 ptHatFilter::ptHatFilter(edm::ParameterSet const& cfg) : token_{consumes(edm::InputTag("genEventScale"))} {
0031 ptHatLowerCut = cfg.getParameter<double>("ptHatLowerCut");
0032 ptHatUpperCut = cfg.getParameter<double>("ptHatUpperCut");
0033 }
0034
0035 void ptHatFilter::beginJob() {
0036 totalEvents = 0;
0037 acceptedEvents = 0;
0038 }
0039
0040
0041 bool ptHatFilter::filter(edm::StreamID, edm::Event& evt, edm::EventSetup const& iSetup) const {
0042 bool result = false;
0043 totalEvents++;
0044 double pt_hat = evt.get(token_);
0045 if (pt_hat > ptHatLowerCut && pt_hat < ptHatUpperCut) {
0046 acceptedEvents++;
0047 result = true;
0048 }
0049 return result;
0050 }
0051
0052 void ptHatFilter::endJob() {
0053 std::cout << "Total Events = " << totalEvents << std::endl;
0054 std::cout << "Accepted Events = " << acceptedEvents << std::endl;
0055 }
0056
0057 #include "FWCore/Framework/interface/MakerMacros.h"
0058 DEFINE_FWK_MODULE(ptHatFilter);