File indexing completed on 2024-04-06 12:13:39
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022 #include <memory>
0023 #include <vector>
0024
0025
0026 #include "FWCore/Framework/interface/global/EDFilter.h"
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/Event.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030
0031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0032 #include "FWCore/Utilities/interface/InputTag.h"
0033
0034 #include "DataFormats/Common/interface/Handle.h"
0035 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0036
0037
0038 class GenHTFilter : public edm::global::EDFilter<> {
0039 public:
0040 explicit GenHTFilter(const edm::ParameterSet&);
0041
0042 private:
0043 bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0044
0045
0046 const edm::EDGetTokenT<reco::GenJetCollection> token_;
0047 const double jetPtCut_, jetEtaCut_, genHTcut_;
0048 };
0049
0050
0051 GenHTFilter::GenHTFilter(const edm::ParameterSet& params)
0052 : token_(consumes<reco::GenJetCollection>(params.getParameter<edm::InputTag>("src"))),
0053 jetPtCut_(params.getParameter<double>("jetPtCut")),
0054 jetEtaCut_(params.getParameter<double>("jetEtaCut")),
0055 genHTcut_(params.getParameter<double>("genHTcut")) {}
0056
0057 bool GenHTFilter::filter(edm::StreamID, edm::Event& evt, const edm::EventSetup& params) const {
0058 using namespace std;
0059 using namespace edm;
0060 using namespace reco;
0061
0062
0063 edm::Handle<reco::GenJetCollection> generatedJets;
0064 evt.getByToken(token_, generatedJets);
0065
0066
0067 double genHT = 0.0;
0068 for (std::vector<reco::GenJet>::const_iterator it = generatedJets->begin(); it != generatedJets->end(); ++it) {
0069 const reco::GenJet& gjet = *it;
0070
0071
0072 if (gjet.pt() > jetPtCut_ && abs(gjet.eta()) < jetEtaCut_) {
0073 genHT += gjet.pt();
0074 }
0075 }
0076 return (genHT > genHTcut_);
0077 }
0078
0079
0080 DEFINE_FWK_MODULE(GenHTFilter);