CMS 3D CMS Logo

SmearedJetProducerT.h
Go to the documentation of this file.
1 #ifndef PAT_SMEAREDJETPRODUCERT_H
2 #define PAT_SMEAREDJETPRODUCERT_H
3 
16 
25 
29 
31 
32 #include <TFile.h>
33 #include <TH1.h>
34 #include <TH1F.h>
35 
36 #include <memory>
37 #include <random>
38 
39 namespace pat {
40  class GenJetMatcher {
41  public:
43  : m_genJetsToken(collector.consumes<reco::GenJetCollection>(cfg.getParameter<edm::InputTag>("genJets"))),
44  m_dR_max(cfg.getParameter<double>("dRMax")),
45  m_dPt_max_factor(cfg.getParameter<double>("dPtMaxFactor")) {
46  // Empty
47  }
48 
50  desc.add<edm::InputTag>("genJets");
51  desc.add<double>("dRMax");
52  desc.add<double>("dPtMaxFactor");
53  }
54 
55  void getTokens(const edm::Event& event) { event.getByToken(m_genJetsToken, m_genJets); }
56 
57  template <class T>
58  const reco::GenJet* match(const T& jet, double resolution) {
60 
61  // Try to find a gen jet matching
62  // dR < m_dR_max
63  // dPt < m_dPt_max_factor * resolution
64 
65  double min_dR = std::numeric_limits<double>::infinity();
66  const reco::GenJet* matched_genJet = nullptr;
67 
68  for (const auto& genJet : genJets) {
69  double dR = deltaR(genJet, jet);
70 
71  if (dR > min_dR)
72  continue;
73 
74  if (dR < m_dR_max) {
75  double dPt = std::abs(genJet.pt() - jet.pt());
76  if (dPt > m_dPt_max_factor * resolution)
77  continue;
78 
79  min_dR = dR;
80  matched_genJet = &genJet;
81  }
82  }
83 
84  return matched_genJet;
85  }
86 
87  private:
90 
91  double m_dR_max;
93  };
94 }; // namespace pat
95 
96 template <typename T>
98  using JetCollection = std::vector<T>;
99 
100 public:
102  : m_enabled(cfg.getParameter<bool>("enabled")),
103  m_useDeterministicSeed(cfg.getParameter<bool>("useDeterministicSeed")),
104  m_debug(cfg.getUntrackedParameter<bool>("debug", false)) {
105  m_jets_token = consumes<JetCollection>(cfg.getParameter<edm::InputTag>("src"));
106 
107  if (m_enabled) {
108  m_rho_token = consumes<double>(cfg.getParameter<edm::InputTag>("rho"));
109 
110  m_use_txt_files = cfg.exists("resolutionFile") && cfg.exists("scaleFactorFile");
111 
112  if (m_use_txt_files) {
113  std::string resolutionFile = cfg.getParameter<edm::FileInPath>("resolutionFile").fullPath();
114  std::string scaleFactorFile = cfg.getParameter<edm::FileInPath>("scaleFactorFile").fullPath();
115 
116  m_resolution_from_file = std::make_unique<JME::JetResolution>(resolutionFile);
117  m_scale_factor_from_file = std::make_unique<JME::JetResolutionScaleFactor>(scaleFactorFile);
118  } else {
119  m_jets_algo = cfg.getParameter<std::string>("algo");
120  m_jets_algo_pt = cfg.getParameter<std::string>("algopt");
123  }
124 
125  std::uint32_t seed = cfg.getParameter<std::uint32_t>("seed");
126  m_random_generator = std::mt19937(seed);
127 
128  bool skipGenMatching = cfg.getParameter<bool>("skipGenMatching");
129  if (!skipGenMatching)
130  m_genJetMatcher = std::make_shared<pat::GenJetMatcher>(cfg, consumesCollector());
131 
132  std::int32_t variation = cfg.getParameter<std::int32_t>("variation");
133  m_uncertaintySource = cfg.getParameter<std::string>("uncertaintySource");
134  m_nomVar = 1;
135  if (variation == 0)
137  else if (variation == 1)
139  else if (variation == -1)
141  else if (variation == 101) {
143  m_nomVar = 1;
144  } else if (variation == -101) {
146  m_nomVar = -1;
147  } else
149  "Invalid value for 'variation' parameter. Only -1, 0, 1 or 101, -101 are supported.");
150  }
151 
152  produces<JetCollection>();
153  }
154 
157 
158  desc.add<edm::InputTag>("src");
159  desc.add<bool>("enabled");
160  desc.add<edm::InputTag>("rho");
161  desc.add<std::int32_t>("variation", 0);
162  desc.add<std::string>("uncertaintySource", "");
163  desc.add<std::uint32_t>("seed", 37428479);
164  desc.add<bool>("skipGenMatching", false);
165  desc.add<bool>("useDeterministicSeed", true);
166  desc.addUntracked<bool>("debug", false);
167 
168  auto source = (edm::ParameterDescription<std::string>("algo", true) and
169  edm::ParameterDescription<std::string>("algopt", true)) xor
170  (edm::ParameterDescription<edm::FileInPath>("resolutionFile", true) and
171  edm::ParameterDescription<edm::FileInPath>("scaleFactorFile", true));
172  desc.addNode(std::move(source));
173 
175 
176  descriptions.addDefault(desc);
177  }
178 
179  void produce(edm::Event& event, const edm::EventSetup& setup) override {
180  edm::Handle<JetCollection> jets_collection;
181  event.getByToken(m_jets_token, jets_collection);
182 
183  // Disable the module when running on real data
184  if (m_enabled && event.isRealData()) {
185  m_enabled = false;
186  m_genJetMatcher.reset();
187  }
188 
190  if (m_enabled)
191  event.getByToken(m_rho_token, rho);
192 
194  JME::JetResolutionScaleFactor resolution_sf;
195 
196  const JetCollection& jets = *jets_collection;
197 
198  if (m_enabled) {
199  if (m_use_txt_files) {
201  resolution_sf = *m_scale_factor_from_file;
202  } else {
205  }
206 
208  unsigned int runNum_uint = static_cast<unsigned int>(event.id().run());
209  unsigned int lumiNum_uint = static_cast<unsigned int>(event.id().luminosityBlock());
210  unsigned int evNum_uint = static_cast<unsigned int>(event.id().event());
211  unsigned int jet0eta = uint32_t(jets.empty() ? 0 : jets[0].eta() / 0.01);
212  std::uint32_t seed = jet0eta + m_nomVar + (lumiNum_uint << 10) + (runNum_uint << 20) + evNum_uint;
213  m_random_generator.seed(seed);
214  }
215  }
216 
217  if (m_genJetMatcher)
218  m_genJetMatcher->getTokens(event);
219 
220  auto smearedJets = std::make_unique<JetCollection>();
221 
222  for (const auto& jet : jets) {
223  if ((!m_enabled) || (jet.pt() == 0)) {
224  // Module disabled or invalid p4. Simply copy the input jet.
225  smearedJets->push_back(jet);
226 
227  continue;
228  }
229 
230  double jet_resolution = resolution.getResolution(
232  double jer_sf = resolution_sf.getScaleFactor({{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}},
235 
236  if (m_debug) {
237  std::cout << "jet: pt: " << jet.pt() << " eta: " << jet.eta() << " phi: " << jet.phi()
238  << " e: " << jet.energy() << std::endl;
239  std::cout << "resolution: " << jet_resolution << std::endl;
240  std::cout << "resolution scale factor: " << jer_sf << std::endl;
241  }
242 
243  const reco::GenJet* genJet = nullptr;
244  if (m_genJetMatcher)
245  genJet = m_genJetMatcher->match(jet, jet.pt() * jet_resolution);
246 
247  double smearFactor = 1.;
248 
249  if (genJet) {
250  /*
251  * Case 1: we have a "good" gen jet matched to the reco jet
252  */
253 
254  if (m_debug) {
255  std::cout << "gen jet: pt: " << genJet->pt() << " eta: " << genJet->eta() << " phi: " << genJet->phi()
256  << " e: " << genJet->energy() << std::endl;
257  }
258 
259  double dPt = jet.pt() - genJet->pt();
260  smearFactor = 1 + m_nomVar * (jer_sf - 1.) * dPt / jet.pt();
261 
262  } else if (jer_sf > 1) {
263  /*
264  * Case 2: we don't have a gen jet. Smear jet pt using a random gaussian variation
265  */
266 
267  double sigma = jet_resolution * std::sqrt(jer_sf * jer_sf - 1);
268  if (m_debug) {
269  std::cout << "gaussian width: " << sigma << std::endl;
270  }
271 
272  std::normal_distribution<> d(0, sigma);
273  smearFactor = 1. + m_nomVar * d(m_random_generator);
274  } else if (m_debug) {
275  std::cout << "Impossible to smear this jet" << std::endl;
276  }
277 
278  if (jet.energy() * smearFactor < MIN_JET_ENERGY) {
279  // Negative or too small smearFactor. We would change direction of the jet
280  // and this is not what we want.
281  // Recompute the smearing factor in order to have jet.energy() == MIN_JET_ENERGY
282  double newSmearFactor = MIN_JET_ENERGY / jet.energy();
283  if (m_debug) {
284  std::cout << "The smearing factor (" << smearFactor << ") is either negative or too small. Fixing it to "
285  << newSmearFactor << " to avoid change of direction." << std::endl;
286  }
287  smearFactor = newSmearFactor;
288  }
289 
290  T smearedJet = jet;
291  smearedJet.scaleEnergy(smearFactor);
292 
293  if (m_debug) {
294  std::cout << "smeared jet (" << smearFactor << "): pt: " << smearedJet.pt() << " eta: " << smearedJet.eta()
295  << " phi: " << smearedJet.phi() << " e: " << smearedJet.energy() << std::endl;
296  }
297 
298  smearedJets->push_back(smearedJet);
299  }
300 
301  // Sort jets by pt
302  std::sort(smearedJets->begin(), smearedJets->end(), jetPtComparator);
303 
304  event.put(std::move(smearedJets));
305  }
306 
307 private:
308  static constexpr const double MIN_JET_ENERGY = 1e-2;
309 
312  bool m_enabled;
320  bool m_debug;
321  std::shared_ptr<pat::GenJetMatcher> m_genJetMatcher;
322 
324  std::unique_ptr<JME::JetResolution> m_resolution_from_file;
325  std::unique_ptr<JME::JetResolutionScaleFactor> m_scale_factor_from_file;
326 
327  std::mt19937 m_random_generator;
328 
330 
331  int m_nomVar;
332 };
333 #endif
edm::EDGetTokenT< JetCollection > m_jets_token
static constexpr const double MIN_JET_ENERGY
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
static const JetResolutionScaleFactor get(const edm::EventSetup &, const Token &)
void getTokens(const edm::Event &event)
double pt() const final
transverse momentum
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::Handle< reco::GenJetCollection > m_genJets
std::vector< GenJet > GenJetCollection
collection of GenJet objects
GenJetMatcher(const edm::ParameterSet &cfg, edm::ConsumesCollector &&collector)
SmearedJetProducerT(const edm::ParameterSet &cfg)
Definition: HeavyIon.h:7
GreaterByPt< T > jetPtComparator
const reco::GenJet * match(const T &jet, double resolution)
JME::JetResolutionScaleFactor::Token m_jets_algo_token
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
float getScaleFactor(const JetParameters &parameters, Variation variation=Variation::NOMINAL, std::string uncertaintySource="") const
T sqrt(T t)
Definition: SSEVec.h:19
const double infinity
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:23
edm::EDGetTokenT< double > m_rho_token
d
Definition: ztail.py:151
std::mt19937 m_random_generator
std::shared_ptr< pat::GenJetMatcher > m_genJetMatcher
fixed size matrix
HLT enums.
std::unique_ptr< JME::JetResolutionScaleFactor > m_scale_factor_from_file
edm::EDGetTokenT< reco::GenJetCollection > m_genJetsToken
long double T
double phi() const final
momentum azimuthal angle
static std::string const source
Definition: EdmProvDump.cc:49
def move(src, dest)
Definition: eostools.py:511
JME::JetResolution::Token m_jets_algo_pt_token
Definition: event.py:1
static const JetResolution get(const edm::EventSetup &, const Token &)
double energy() const final
energy
double eta() const final
momentum pseudorapidity