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.reset(new JME::JetResolution(resolutionFile));
117  m_scale_factor_from_file.reset(new JME::JetResolutionScaleFactor(scaleFactorFile));
118  } else {
119  m_jets_algo = cfg.getParameter<std::string>("algo");
120  m_jets_algo_pt = cfg.getParameter<std::string>("algopt");
121  }
122 
123  std::uint32_t seed = cfg.getParameter<std::uint32_t>("seed");
124  m_random_generator = std::mt19937(seed);
125 
126  bool skipGenMatching = cfg.getParameter<bool>("skipGenMatching");
127  if (!skipGenMatching)
128  m_genJetMatcher = std::make_shared<pat::GenJetMatcher>(cfg, consumesCollector());
129 
130  std::int32_t variation = cfg.getParameter<std::int32_t>("variation");
131  m_nomVar = 1;
132  if (variation == 0)
133  m_systematic_variation = Variation::NOMINAL;
134  else if (variation == 1)
135  m_systematic_variation = Variation::UP;
136  else if (variation == -1)
137  m_systematic_variation = Variation::DOWN;
138  else if (variation == 101) {
139  m_systematic_variation = Variation::NOMINAL;
140  m_nomVar = 1;
141  } else if (variation == -101) {
142  m_systematic_variation = Variation::NOMINAL;
143  m_nomVar = -1;
144  } else
146  "Invalid value for 'variation' parameter. Only -1, 0, 1 or 101, -101 are supported.");
147  }
148 
149  produces<JetCollection>();
150  }
151 
154 
155  desc.add<edm::InputTag>("src");
156  desc.add<bool>("enabled");
157  desc.add<edm::InputTag>("rho");
158  desc.add<std::int32_t>("variation", 0);
159  desc.add<std::uint32_t>("seed", 37428479);
160  desc.add<bool>("skipGenMatching", false);
161  desc.add<bool>("useDeterministicSeed", true);
162  desc.addUntracked<bool>("debug", false);
163 
164  auto source = (edm::ParameterDescription<std::string>("algo", true) and
165  edm::ParameterDescription<std::string>("algopt", true)) xor
166  (edm::ParameterDescription<edm::FileInPath>("resolutionFile", true) and
167  edm::ParameterDescription<edm::FileInPath>("scaleFactorFile", true));
168  desc.addNode(std::move(source));
169 
171 
172  descriptions.addDefault(desc);
173  }
174 
175  void produce(edm::Event& event, const edm::EventSetup& setup) override {
176  edm::Handle<JetCollection> jets_collection;
177  event.getByToken(m_jets_token, jets_collection);
178 
179  // Disable the module when running on real data
180  if (m_enabled && event.isRealData()) {
181  m_enabled = false;
182  m_genJetMatcher.reset();
183  }
184 
186  if (m_enabled)
187  event.getByToken(m_rho_token, rho);
188 
190  JME::JetResolutionScaleFactor resolution_sf;
191 
192  const JetCollection& jets = *jets_collection;
193 
194  if (m_enabled) {
195  if (m_use_txt_files) {
196  resolution = *m_resolution_from_file;
197  resolution_sf = *m_scale_factor_from_file;
198  } else {
199  resolution = JME::JetResolution::get(setup, m_jets_algo_pt);
200  resolution_sf = JME::JetResolutionScaleFactor::get(setup, m_jets_algo);
201  }
202 
203  if (m_useDeterministicSeed) {
204  unsigned int runNum_uint = static_cast<unsigned int>(event.id().run());
205  unsigned int lumiNum_uint = static_cast<unsigned int>(event.id().luminosityBlock());
206  unsigned int evNum_uint = static_cast<unsigned int>(event.id().event());
207  unsigned int jet0eta = uint32_t(jets.empty() ? 0 : jets[0].eta() / 0.01);
208  std::uint32_t seed = jet0eta + m_nomVar + (lumiNum_uint << 10) + (runNum_uint << 20) + evNum_uint;
209  m_random_generator.seed(seed);
210  }
211  }
212 
213  if (m_genJetMatcher)
214  m_genJetMatcher->getTokens(event);
215 
216  auto smearedJets = std::make_unique<JetCollection>();
217 
218  for (const auto& jet : jets) {
219  if ((!m_enabled) || (jet.pt() == 0)) {
220  // Module disabled or invalid p4. Simply copy the input jet.
221  smearedJets->push_back(jet);
222 
223  continue;
224  }
225 
226  double jet_resolution = resolution.getResolution(
227  {{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}, {JME::Binning::Rho, *rho}});
228  double jer_sf = resolution_sf.getScaleFactor({{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}},
229  m_systematic_variation);
230 
231  if (m_debug) {
232  std::cout << "jet: pt: " << jet.pt() << " eta: " << jet.eta() << " phi: " << jet.phi()
233  << " e: " << jet.energy() << std::endl;
234  std::cout << "resolution: " << jet_resolution << std::endl;
235  std::cout << "resolution scale factor: " << jer_sf << std::endl;
236  }
237 
238  const reco::GenJet* genJet = nullptr;
239  if (m_genJetMatcher)
240  genJet = m_genJetMatcher->match(jet, jet.pt() * jet_resolution);
241 
242  double smearFactor = 1.;
243 
244  if (genJet) {
245  /*
246  * Case 1: we have a "good" gen jet matched to the reco jet
247  */
248 
249  if (m_debug) {
250  std::cout << "gen jet: pt: " << genJet->pt() << " eta: " << genJet->eta() << " phi: " << genJet->phi()
251  << " e: " << genJet->energy() << std::endl;
252  }
253 
254  double dPt = jet.pt() - genJet->pt();
255  smearFactor = 1 + m_nomVar * (jer_sf - 1.) * dPt / jet.pt();
256 
257  } else if (jer_sf > 1) {
258  /*
259  * Case 2: we don't have a gen jet. Smear jet pt using a random gaussian variation
260  */
261 
262  double sigma = jet_resolution * std::sqrt(jer_sf * jer_sf - 1);
263  if (m_debug) {
264  std::cout << "gaussian width: " << sigma << std::endl;
265  }
266 
267  std::normal_distribution<> d(0, sigma);
268  smearFactor = 1. + m_nomVar * d(m_random_generator);
269  } else if (m_debug) {
270  std::cout << "Impossible to smear this jet" << std::endl;
271  }
272 
273  if (jet.energy() * smearFactor < MIN_JET_ENERGY) {
274  // Negative or too small smearFactor. We would change direction of the jet
275  // and this is not what we want.
276  // Recompute the smearing factor in order to have jet.energy() == MIN_JET_ENERGY
277  double newSmearFactor = MIN_JET_ENERGY / jet.energy();
278  if (m_debug) {
279  std::cout << "The smearing factor (" << smearFactor << ") is either negative or too small. Fixing it to "
280  << newSmearFactor << " to avoid change of direction." << std::endl;
281  }
282  smearFactor = newSmearFactor;
283  }
284 
285  T smearedJet = jet;
286  smearedJet.scaleEnergy(smearFactor);
287 
288  if (m_debug) {
289  std::cout << "smeared jet (" << smearFactor << "): pt: " << smearedJet.pt() << " eta: " << smearedJet.eta()
290  << " phi: " << smearedJet.phi() << " e: " << smearedJet.energy() << std::endl;
291  }
292 
293  smearedJets->push_back(smearedJet);
294  }
295 
296  // Sort jets by pt
297  std::sort(smearedJets->begin(), smearedJets->end(), jetPtComparator);
298 
299  event.put(std::move(smearedJets));
300  }
301 
302 private:
303  static constexpr const double MIN_JET_ENERGY = 1e-2;
304 
307  bool m_enabled;
312  bool m_debug;
313  std::shared_ptr<pat::GenJetMatcher> m_genJetMatcher;
314 
316  std::unique_ptr<JME::JetResolution> m_resolution_from_file;
317  std::unique_ptr<JME::JetResolutionScaleFactor> m_scale_factor_from_file;
318 
319  std::mt19937 m_random_generator;
320 
322 
323  int m_nomVar;
324 };
325 #endif
edm::EDGetTokenT< JetCollection > m_jets_token
T getParameter(std::string const &) const
float getResolution(const JetParameters &parameters) 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 &parameterName) const
checks if a parameter exists
SmearedJetProducerT(const edm::ParameterSet &cfg)
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
double pt() const final
transverse momentum
bool isRealData() const
Definition: EventBase.h:62
Definition: HeavyIon.h:7
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
T sqrt(T t)
Definition: SSEVec.h:19
double energy() const final
energy
const double infinity
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:23
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< double > m_rho_token
d
Definition: ztail.py:151
std::mt19937 m_random_generator
float getScaleFactor(const JetParameters &parameters, Variation variation=Variation::NOMINAL) const
static const JetResolutionScaleFactor get(const edm::EventSetup &, const std::string &)
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:47
def move(src, dest)
Definition: eostools.py:511
#define constexpr
Definition: event.py:1