Skip to content

Commit

Permalink
Merge pull request #27794 from ohindrichs/toplhefilter
Browse files Browse the repository at this point in the history
Add TopLHEFilter
  • Loading branch information
cmsbuild committed Aug 26, 2019
2 parents 6c7535c + 1432426 commit 0920e9c
Show file tree
Hide file tree
Showing 3 changed files with 75 additions and 0 deletions.
1 change: 1 addition & 0 deletions GeneratorInterface/GenFilters/BuildFile.xml
@@ -1,3 +1,4 @@
<use name="fastjet"/>
<use name="boost"/>
<use name="FWCore/PluginManager"/>
<use name="FWCore/ParameterSet"/>
Expand Down
8 changes: 8 additions & 0 deletions GeneratorInterface/GenFilters/python/LHEJetFilter_cfi.py
@@ -0,0 +1,8 @@
import FWCore.ParameterSet.Config as cms

LHEJetFilter = cms.EDFilter('LHEJetFilter',
jetPtMin = cms.double(350.),
jetR = cms.double(0.8),
src = cms.InputTag('externalLHEProducer')
)

66 changes: 66 additions & 0 deletions GeneratorInterface/GenFilters/src/LHEJetFilter.cc
@@ -0,0 +1,66 @@
#include <vector>

#include "FWCore/Framework/interface/global/EDFilter.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/InputTag.h"

#include "DataFormats/Common/interface/Handle.h"

#include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"

#include "fastjet/PseudoJet.hh"
#include "fastjet/JetDefinition.hh"
#include "fastjet/ClusterSequence.hh"
#include "fastjet/Selector.hh"

using namespace std;
using namespace fastjet;

class LHEJetFilter : public edm::global::EDFilter<> {
public:
explicit LHEJetFilter(const edm::ParameterSet&);
~LHEJetFilter() override {}

private:
bool filter(edm::StreamID strid, edm::Event& evt, const edm::EventSetup& params) const override;

edm::EDGetTokenT<LHEEventProduct> tokenLHEEvent_;
double jetPtMin_;
JetDefinition jetdef_;
};

LHEJetFilter::LHEJetFilter(const edm::ParameterSet& params)
: tokenLHEEvent_(consumes<LHEEventProduct>(params.getParameter<edm::InputTag>("src"))),
jetPtMin_(params.getParameter<double>("jetPtMin")),
jetdef_(antikt_algorithm, params.getParameter<double>("jetR")) {}

bool LHEJetFilter::filter(edm::StreamID strid, edm::Event& evt, const edm::EventSetup& params) const {
edm::Handle<LHEEventProduct> lheinfo;
evt.getByToken(tokenLHEEvent_, lheinfo);

if (!lheinfo.isValid()) {
return true;
}

vector<PseudoJet> jetconsts;
jetconsts.reserve(10);
const lhef::HEPEUP& hepeup = lheinfo->hepeup();
for (size_t p = 0; p < hepeup.IDUP.size(); ++p) {
if (hepeup.ISTUP[p] == 1) {
const lhef::HEPEUP::FiveVector& mom = hepeup.PUP[p];
jetconsts.emplace_back(mom[0], mom[1], mom[2], mom[3]);
}
}

ClusterSequence cs(jetconsts, jetdef_);
vector<PseudoJet> jets = cs.inclusive_jets(jetPtMin_);

return !jets.empty();
}

// Define module as a plug-in
DEFINE_FWK_MODULE(LHEJetFilter);

0 comments on commit 0920e9c

Please sign in to comment.