1 #ifndef PAT_SMEAREDJETPRODUCERT_H 2 #define PAT_SMEAREDJETPRODUCERT_H 44 m_dR_max(cfg.getParameter<double>(
"dRMax")),
51 desc.
add<
double>(
"dRMax");
52 desc.
add<
double>(
"dPtMaxFactor");
70 for (
const auto& genJet: genJets) {
77 double dPt =
std::abs(genJet.pt() - jet.pt());
82 matched_genJet = &genJet;
86 return matched_genJet;
105 m_enabled(cfg.getParameter<bool>(
"enabled")),
106 m_useDeterministicSeed(cfg.getParameter<bool>(
"useDeterministicSeed")),
107 m_debug(cfg.getUntrackedParameter<bool>(
"debug",
false)) {
114 m_use_txt_files = cfg.
exists(
"resolutionFile") && cfg.
exists(
"scaleFactorFile");
116 if (m_use_txt_files) {
128 m_random_generator = std::mt19937(seed);
131 if (! skipGenMatching)
132 m_genJetMatcher = std::make_shared<pat::GenJetMatcher>(
cfg, consumesCollector());
134 std::int32_t
variation = cfg.getParameter<std::int32_t>(
"variation");
138 else if (variation == 1)
140 else if (variation == -1)
142 else if (variation == 101) {
146 else if (variation == -101) {
154 produces<JetCollection>();
161 desc.
add<
bool>(
"enabled");
163 desc.
add<std::int32_t>(
"variation", 0);
164 desc.
add<std::uint32_t>(
"seed", 37428479);
165 desc.
add<
bool>(
"skipGenMatching",
false);
166 desc.
add<
bool>(
"useDeterministicSeed",
true);
182 event.getByToken(m_jets_token, jets_collection);
187 m_genJetMatcher.reset();
192 event.getByToken(m_rho_token, rho);
200 if (m_use_txt_files) {
201 resolution = *m_resolution_from_file;
202 resolution_sf = *m_scale_factor_from_file;
208 if(m_useDeterministicSeed) {
209 unsigned int runNum_uint = static_cast <
unsigned int> (
event.id().run());
210 unsigned int lumiNum_uint = static_cast <
unsigned int> (
event.id().luminosityBlock());
211 unsigned int evNum_uint = static_cast <
unsigned int> (
event.id().event());
212 unsigned int jet0eta = uint32_t(jets.empty() ? 0 : jets[0].eta()/0.01);
213 std::uint32_t
seed = jet0eta + m_nomVar + (lumiNum_uint<<10) + (runNum_uint<<20) + evNum_uint;
214 m_random_generator.seed(seed);
219 m_genJetMatcher->getTokens(event);
221 auto smearedJets = std::make_unique<JetCollection>();
223 for (
const auto&
jet: jets) {
225 if ((! m_enabled) || (
jet.pt() == 0)) {
227 smearedJets->push_back(
jet);
236 std::cout <<
"jet: pt: " <<
jet.pt() <<
" eta: " <<
jet.eta() <<
" phi: " <<
jet.phi() <<
" e: " <<
jet.energy() << std::endl;
237 std::cout <<
"resolution: " << jet_resolution << std::endl;
238 std::cout <<
"resolution scale factor: " << jer_sf << std::endl;
243 genJet = m_genJetMatcher->match(
jet,
jet.pt() * jet_resolution);
245 double smearFactor = 1.;
253 std::cout <<
"gen jet: pt: " << genJet->
pt() <<
" eta: " << genJet->
eta() <<
" phi: " << genJet->
phi() <<
" e: " << genJet->
energy() << std::endl;
256 double dPt =
jet.pt() - genJet->
pt();
257 smearFactor = 1 + m_nomVar*(jer_sf - 1.) * dPt /
jet.pt();
259 }
else if (jer_sf > 1) {
264 double sigma = jet_resolution *
std::sqrt(jer_sf * jer_sf - 1);
266 std::cout <<
"gaussian width: " << sigma << std::endl;
269 std::normal_distribution<>
d(0, sigma);
270 smearFactor = 1. + m_nomVar*
d(m_random_generator);
271 }
else if (m_debug) {
272 std::cout <<
"Impossible to smear this jet" << std::endl;
275 if (
jet.energy() * smearFactor < MIN_JET_ENERGY) {
279 double newSmearFactor = MIN_JET_ENERGY /
jet.energy();
281 std::cout <<
"The smearing factor (" << smearFactor <<
") is either negative or too small. Fixing it to " << newSmearFactor <<
" to avoid change of direction." << std::endl;
283 smearFactor = newSmearFactor;
287 smearedJet.scaleEnergy(smearFactor);
290 std::cout <<
"smeared jet (" << smearFactor <<
"): pt: " << smearedJet.pt() <<
" eta: " << smearedJet.eta() <<
" phi: " << smearedJet.phi() <<
" e: " << smearedJet.energy() << std::endl;
293 smearedJets->push_back(smearedJet);
297 std::sort(smearedJets->begin(), smearedJets->end(), jetPtComparator);
edm::EDGetTokenT< JetCollection > m_jets_token
T getParameter(std::string const &) const
virtual double pt() const final
transverse momentum
float getResolution(const JetParameters ¶meters) const
void getTokens(const edm::Event &event)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static const JetResolution get(const edm::EventSetup &, const std::string &)
virtual double eta() const final
momentum pseudorapidity
edm::Handle< reco::GenJetCollection > m_genJets
std::vector< GenJet > GenJetCollection
collection of GenJet objects
virtual void produce(edm::Event &event, const edm::EventSetup &setup) override
GenJetMatcher(const edm::ParameterSet &cfg, edm::ConsumesCollector &&collector)
bool exists(std::string const ¶meterName) const
checks if a parameter exists
SmearedJetProducerT(const edm::ParameterSet &cfg)
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
def setup(process, global_tag, zero_tesla=False)
GreaterByPt< T > jetPtComparator
const reco::GenJet * match(const T &jet, double resolution)
virtual double phi() const final
momentum azimuthal angle
std::unique_ptr< JME::JetResolution > m_resolution_from_file
static void fillDescriptions(edm::ParameterSetDescription &desc)
void addDefault(ParameterSetDescription const &psetDescription)
std::vector< T > JetCollection
std::string m_jets_algo_pt
bool m_useDeterministicSeed
virtual double energy() const final
energy
Abs< T >::type abs(const T &t)
Jets made from MC generator particles.
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< double > m_rho_token
double deltaR(double eta1, double eta2, double phi1, double phi2)
std::mt19937 m_random_generator
float getScaleFactor(const JetParameters ¶meters, Variation variation=Variation::NOMINAL) const
Variation m_systematic_variation
static const JetResolutionScaleFactor get(const edm::EventSetup &, const std::string &)
std::shared_ptr< pat::GenJetMatcher > m_genJetMatcher
std::unique_ptr< JME::JetResolutionScaleFactor > m_scale_factor_from_file
edm::EDGetTokenT< reco::GenJetCollection > m_genJetsToken
static std::string const source