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");
68 for (
const auto& genJet : genJets) {
75 double dPt =
std::abs(genJet.pt() - jet.pt());
80 matched_genJet = &genJet;
84 return matched_genJet;
102 : m_enabled(cfg.getParameter<
bool>(
"enabled")),
103 m_useDeterministicSeed(cfg.getParameter<
bool>(
"useDeterministicSeed")),
104 m_debug(cfg.getUntrackedParameter<
bool>(
"debug",
false)) {
110 m_use_txt_files = cfg.
exists(
"resolutionFile") && cfg.
exists(
"scaleFactorFile");
112 if (m_use_txt_files) {
124 m_random_generator = std::mt19937(seed);
127 if (!skipGenMatching)
128 m_genJetMatcher = std::make_shared<pat::GenJetMatcher>(
cfg, consumesCollector());
130 std::int32_t
variation = cfg.getParameter<std::int32_t>(
"variation");
131 m_uncertaintySource = cfg.getParameter<
std::string>(
"uncertaintySource");
135 else if (variation == 1)
137 else if (variation == -1)
139 else if (variation == 101) {
142 }
else if (variation == -101) {
147 "Invalid value for 'variation' parameter. Only -1, 0, 1 or 101, -101 are supported.");
150 produces<JetCollection>();
157 desc.
add<
bool>(
"enabled");
159 desc.
add<std::int32_t>(
"variation", 0);
161 desc.
add<std::uint32_t>(
"seed", 37428479);
162 desc.
add<
bool>(
"skipGenMatching",
false);
163 desc.
add<
bool>(
"useDeterministicSeed",
true);
179 event.getByToken(m_jets_token, jets_collection);
184 m_genJetMatcher.reset();
189 event.getByToken(m_rho_token, rho);
197 if (m_use_txt_files) {
198 resolution = *m_resolution_from_file;
199 resolution_sf = *m_scale_factor_from_file;
205 if (m_useDeterministicSeed) {
206 unsigned int runNum_uint =
static_cast<unsigned int>(
event.id().run());
207 unsigned int lumiNum_uint =
static_cast<unsigned int>(
event.id().luminosityBlock());
208 unsigned int evNum_uint =
static_cast<unsigned int>(
event.id().event());
209 unsigned int jet0eta = uint32_t(jets.empty() ? 0 : jets[0].eta() / 0.01);
210 std::uint32_t
seed = jet0eta + m_nomVar + (lumiNum_uint << 10) + (runNum_uint << 20) + evNum_uint;
211 m_random_generator.seed(seed);
216 m_genJetMatcher->getTokens(event);
218 auto smearedJets = std::make_unique<JetCollection>();
220 for (
const auto&
jet : jets) {
221 if ((!m_enabled) || (
jet.pt() == 0)) {
223 smearedJets->push_back(
jet);
231 m_systematic_variation,
232 m_uncertaintySource);
234 std::cout <<
"jet: pt: " <<
jet.pt() <<
" eta: " <<
jet.eta() <<
" phi: " <<
jet.phi()
235 <<
" e: " <<
jet.energy() << std::endl;
236 std::cout <<
"resolution: " << jet_resolution << std::endl;
237 std::cout <<
"resolution scale factor: " << jer_sf << std::endl;
242 genJet = m_genJetMatcher->match(
jet,
jet.pt() * jet_resolution);
244 double smearFactor = 1.;
252 std::cout <<
"gen jet: pt: " << genJet->
pt() <<
" eta: " << genJet->
eta() <<
" phi: " << genJet->
phi()
253 <<
" 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 " 282 << newSmearFactor <<
" to avoid change of direction." << std::endl;
284 smearFactor = newSmearFactor;
288 smearedJet.scaleEnergy(smearFactor);
291 std::cout <<
"smeared jet (" << smearFactor <<
"): pt: " << smearedJet.pt() <<
" eta: " << smearedJet.eta()
292 <<
" phi: " << smearedJet.phi() <<
" e: " << smearedJet.energy() << std::endl;
295 smearedJets->push_back(smearedJet);
299 std::sort(smearedJets->begin(), smearedJets->end(), jetPtComparator);
edm::EDGetTokenT< JetCollection > m_jets_token
T getParameter(std::string const &) const
float getResolution(const JetParameters ¶meters) const
double eta() const final
momentum pseudorapidity
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 &)
edm::Handle< reco::GenJetCollection > m_genJets
std::vector< GenJet > GenJetCollection
collection of GenJet objects
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)
double pt() const final
transverse momentum
std::string m_uncertaintySource
GreaterByPt< T > jetPtComparator
const reco::GenJet * match(const T &jet, double resolution)
std::unique_ptr< JME::JetResolution > m_resolution_from_file
void produce(edm::Event &event, const edm::EventSetup &setup) override
static void fillDescriptions(edm::ParameterSetDescription &desc)
void addDefault(ParameterSetDescription const &psetDescription)
std::vector< T > JetCollection
std::string m_jets_algo_pt
bool m_useDeterministicSeed
double energy() const final
energy
Abs< T >::type abs(const T &t)
Jets made from MC generator particles.
float getScaleFactor(const JetParameters ¶meters, Variation variation=Variation::NOMINAL, std::string uncertaintySource="") const
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
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
double phi() const final
momentum azimuthal angle
static std::string const source