File indexing completed on 2024-04-06 12:23:25
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 <string>
0024
0025
0026 #include "FWCore/Framework/interface/Frameworkfwd.h"
0027 #include "FWCore/Framework/interface/global/EDFilter.h"
0028
0029 #include "FWCore/Framework/interface/Event.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033
0034 #include "DataFormats/JetReco/interface/GenJet.h"
0035
0036 #include "PhysicsTools/JetMCUtils/interface/JetMCTag.h"
0037
0038 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0039
0040
0041
0042
0043
0044 class HFFilter : public edm::global::EDFilter<> {
0045 public:
0046 explicit HFFilter(const edm::ParameterSet&);
0047
0048 bool filter(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0049
0050 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0051
0052 private:
0053
0054 edm::EDGetTokenT<std::vector<reco::GenJet>> genJetsCollToken_;
0055 double ptMin_;
0056 double etaMax_;
0057 };
0058 using namespace std;
0059
0060
0061
0062
0063 HFFilter::HFFilter(const edm::ParameterSet& iConfig) {
0064 genJetsCollToken_ = consumes<vector<reco::GenJet>>(iConfig.getParameter<edm::InputTag>("genJetsCollName"));
0065 ptMin_ = iConfig.getParameter<double>("ptMin");
0066 etaMax_ = iConfig.getParameter<double>("etaMax");
0067 }
0068
0069 void HFFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0070 edm::ParameterSetDescription desc;
0071 desc.add<edm::InputTag>("genJetsCollName");
0072 desc.add<double>("ptMin");
0073 desc.add<double>("etaMax");
0074
0075 descriptions.addDefault(desc);
0076 }
0077
0078
0079
0080
0081
0082
0083
0084 bool HFFilter::filter(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0085
0086 using namespace edm;
0087 using namespace reco;
0088 Handle<std::vector<GenJet>> hGenJets;
0089 iEvent.getByToken(genJetsCollToken_, hGenJets);
0090
0091
0092 vector<GenJet>::const_iterator ijet = hGenJets->begin();
0093 vector<GenJet>::const_iterator end = hGenJets->end();
0094 for (; ijet != end; ++ijet) {
0095
0096 if (ijet->pt() < ptMin_ || fabs(ijet->eta()) > etaMax_)
0097 continue;
0098
0099
0100 vector<const GenParticle*> particles = ijet->getGenConstituents();
0101
0102
0103 vector<const GenParticle*>::const_iterator genit = particles.begin();
0104 vector<const GenParticle*>::const_iterator genend = particles.end();
0105 for (; genit != genend; ++genit) {
0106
0107 const GenParticle& genitref = **genit;
0108 if (JetMCTagUtils::decayFromBHadron(genitref) || JetMCTagUtils::decayFromCHadron(genitref)) {
0109 return true;
0110 }
0111 }
0112 }
0113
0114 return false;
0115 }
0116
0117
0118 DEFINE_FWK_MODULE(HFFilter);