CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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) {
56  event.getByToken(m_genJetsToken, m_genJets);
57  }
58 
59  template<class T>
60  const reco::GenJet* match(const T& jet, double resolution) {
61  const reco::GenJetCollection& genJets = *m_genJets;
62 
63  // Try to find a gen jet matching
64  // dR < m_dR_max
65  // dPt < m_dPt_max_factor * resolution
66 
67  double min_dR = std::numeric_limits<double>::infinity();
68  const reco::GenJet* matched_genJet = nullptr;
69 
70  for (const auto& genJet: genJets) {
71  double dR = deltaR(genJet, jet);
72 
73  if (dR > min_dR)
74  continue;
75 
76  if (dR < m_dR_max) {
77  double dPt = std::abs(genJet.pt() - jet.pt());
78  if (dPt > m_dPt_max_factor * resolution)
79  continue;
80 
81  min_dR = dR;
82  matched_genJet = &genJet;
83  }
84  }
85 
86  return matched_genJet;
87  }
88 
89  private:
92 
93  double m_dR_max;
95  };
96 };
97 
98 template <typename T>
100 
101  using JetCollection = std::vector<T>;
102 
103  public:
105  m_enabled(cfg.getParameter<bool>("enabled")),
106  m_useDeterministicSeed(cfg.getParameter<bool>("useDeterministicSeed")),
107  m_debug(cfg.getUntrackedParameter<bool>("debug", false)) {
108 
109  m_jets_token = consumes<JetCollection>(cfg.getParameter<edm::InputTag>("src"));
110 
111  if (m_enabled) {
112  m_rho_token = consumes<double>(cfg.getParameter<edm::InputTag>("rho"));
113 
114  m_use_txt_files = cfg.exists("resolutionFile") && cfg.exists("scaleFactorFile");
115 
116  if (m_use_txt_files) {
117  std::string resolutionFile = cfg.getParameter<edm::FileInPath>("resolutionFile").fullPath();
118  std::string scaleFactorFile = cfg.getParameter<edm::FileInPath>("scaleFactorFile").fullPath();
119 
120  m_resolution_from_file.reset(new JME::JetResolution(resolutionFile));
121  m_scale_factor_from_file.reset(new JME::JetResolutionScaleFactor(scaleFactorFile));
122  } else {
123  m_jets_algo = cfg.getParameter<std::string>("algo");
124  m_jets_algo_pt = cfg.getParameter<std::string>("algopt");
125  }
126 
127  std::uint32_t seed = cfg.getParameter<std::uint32_t>("seed");
128  m_random_generator = std::mt19937(seed);
129 
130  bool skipGenMatching = cfg.getParameter<bool>("skipGenMatching");
131  if (! skipGenMatching)
132  m_genJetMatcher = std::make_shared<pat::GenJetMatcher>(cfg, consumesCollector());
133 
134  std::int32_t variation = cfg.getParameter<std::int32_t>("variation");
135  if (variation == 0)
137  else if (variation == 1)
139  else if (variation == -1)
141  else
142  throw edm::Exception(edm::errors::ConfigFileReadError, "Invalid value for 'variation' parameter. Only -1, 0 or 1 are supported.");
143  }
144 
145  produces<JetCollection>();
146  }
147 
150 
151  desc.add<edm::InputTag>("src");
152  desc.add<bool>("enabled");
153  desc.add<edm::InputTag>("rho");
154  desc.add<std::int32_t>("variation", 0);
155  desc.add<std::uint32_t>("seed", 37428479);
156  desc.add<bool>("skipGenMatching", false);
157  desc.add<bool>("useDeterministicSeed", true);
158  desc.addUntracked<bool>("debug", false);
159 
160  auto source =
162  (edm::ParameterDescription<edm::FileInPath>("resolutionFile", true) and edm::ParameterDescription<edm::FileInPath>("scaleFactorFile", true));
163  desc.addNode(source);
164 
166 
167  descriptions.addDefault(desc);
168  }
169 
170  virtual void produce(edm::Event& event, const edm::EventSetup& setup) override {
171 
172  edm::Handle<JetCollection> jets_collection;
173  event.getByToken(m_jets_token, jets_collection);
174 
175  // Disable the module when running on real data
176  if (m_enabled && event.isRealData()) {
177  m_enabled = false;
178  m_genJetMatcher.reset();
179  }
180 
182  if (m_enabled)
183  event.getByToken(m_rho_token, rho);
184 
186  JME::JetResolutionScaleFactor resolution_sf;
187 
188  const JetCollection& jets = *jets_collection;
189 
190  if (m_enabled) {
191  if (m_use_txt_files) {
192  resolution = *m_resolution_from_file;
193  resolution_sf = *m_scale_factor_from_file;
194  } else {
195  resolution = JME::JetResolution::get(setup, m_jets_algo_pt);
196  resolution_sf = JME::JetResolutionScaleFactor::get(setup, m_jets_algo);
197  }
198 
200  unsigned int runNum_uint = static_cast <unsigned int> (event.id().run());
201  unsigned int lumiNum_uint = static_cast <unsigned int> (event.id().luminosityBlock());
202  unsigned int evNum_uint = static_cast <unsigned int> (event.id().event());
203  unsigned int jet0eta = uint32_t(jets.empty() ? 0 : jets[0].eta()/0.01);
204  std::uint32_t seed = jet0eta + (lumiNum_uint<<10) + (runNum_uint<<20) + evNum_uint;
205  m_random_generator.seed(seed);
206  }
207  }
208 
209  if (m_genJetMatcher)
210  m_genJetMatcher->getTokens(event);
211 
212  std::unique_ptr<JetCollection> smearedJets(new JetCollection());
213 
214  for (const auto& jet: jets) {
215 
216  if ((! m_enabled) || (jet.pt() == 0)) {
217  // Module disabled or invalid p4. Simply copy the input jet.
218  smearedJets->push_back(jet);
219 
220  continue;
221  }
222 
223  double jet_resolution = resolution.getResolution({{JME::Binning::JetPt, jet.pt()}, {JME::Binning::JetEta, jet.eta()}, {JME::Binning::Rho, *rho}});
224  double jer_sf = resolution_sf.getScaleFactor({{JME::Binning::JetEta, jet.eta()}}, m_systematic_variation);
225 
226  if (m_debug) {
227  std::cout << "jet: pt: " << jet.pt() << " eta: " << jet.eta() << " phi: " << jet.phi() << " e: " << jet.energy() << std::endl;
228  std::cout << "resolution: " << jet_resolution << std::endl;
229  std::cout << "resolution scale factor: " << jer_sf << std::endl;
230  }
231 
232  const reco::GenJet* genJet = nullptr;
233  if (m_genJetMatcher)
234  genJet = m_genJetMatcher->match(jet, jet.pt() * jet_resolution);
235 
236  double smearFactor = 1.;
237 
238  if (genJet) {
239  /*
240  * Case 1: we have a "good" gen jet matched to the reco jet
241  */
242 
243  if (m_debug) {
244  std::cout << "gen jet: pt: " << genJet->pt() << " eta: " << genJet->eta() << " phi: " << genJet->phi() << " e: " << genJet->energy() << std::endl;
245  }
246 
247  double dPt = jet.pt() - genJet->pt();
248  smearFactor = 1 + (jer_sf - 1.) * dPt / jet.pt();
249 
250  } else if (jer_sf > 1) {
251  /*
252  * Case 2: we don't have a gen jet. Smear jet pt using a random gaussian variation
253  */
254 
255  double sigma = jet_resolution * std::sqrt(jer_sf * jer_sf - 1);
256  if (m_debug) {
257  std::cout << "gaussian width: " << sigma << std::endl;
258  }
259 
260  std::normal_distribution<> d(0, sigma);
261  smearFactor = 1. + d(m_random_generator);
262  } else if (m_debug) {
263  std::cout << "Impossible to smear this jet" << std::endl;
264  }
265 
266  if (jet.energy() * smearFactor < MIN_JET_ENERGY) {
267  // Negative or too small smearFactor. We would change direction of the jet
268  // and this is not what we want.
269  // Recompute the smearing factor in order to have jet.energy() == MIN_JET_ENERGY
270  double newSmearFactor = MIN_JET_ENERGY / jet.energy();
271  if (m_debug) {
272  std::cout << "The smearing factor (" << smearFactor << ") is either negative or too small. Fixing it to " << newSmearFactor << " to avoid change of direction." << std::endl;
273  }
274  smearFactor = newSmearFactor;
275  }
276 
277  T smearedJet = jet;
278  smearedJet.scaleEnergy(smearFactor);
279 
280  if (m_debug) {
281  std::cout << "smeared jet (" << smearFactor << "): pt: " << smearedJet.pt() << " eta: " << smearedJet.eta() << " phi: " << smearedJet.phi() << " e: " << smearedJet.energy() << std::endl;
282  }
283 
284  smearedJets->push_back(smearedJet);
285  }
286 
287  // Sort jets by pt
288  std::sort(smearedJets->begin(), smearedJets->end(), jetPtComparator);
289 
290  event.put(std::move(smearedJets));
291  }
292 
293  private:
294  static constexpr const double MIN_JET_ENERGY = 1e-2;
295 
298  bool m_enabled;
303  bool m_debug;
304  std::shared_ptr<pat::GenJetMatcher> m_genJetMatcher;
305 
307  std::unique_ptr<JME::JetResolution> m_resolution_from_file;
308  std::unique_ptr<JME::JetResolutionScaleFactor> m_scale_factor_from_file;
309 
310  std::mt19937 m_random_generator;
311 
313 };
314 #endif
edm::EDGetTokenT< JetCollection > m_jets_token
T getParameter(std::string const &) const
float getResolution(const JetParameters &parameters) const
tuple cfg
Definition: looper.py:293
void getTokens(const edm::Event &event)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual double energy() const final
energy
static const JetResolution get(const edm::EventSetup &, const std::string &)
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 &parameterName) const
checks if a parameter exists
SmearedJetProducerT(const edm::ParameterSet &cfg)
virtual double phi() const final
momentum azimuthal angle
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
bool isRealData() const
Definition: EventBase.h:63
#define constexpr
GreaterByPt< T > jetPtComparator
const reco::GenJet * match(const T &jet, double resolution)
std::unique_ptr< JME::JetResolution > m_resolution_from_file
tuple d
Definition: ztail.py:151
static void fillDescriptions(edm::ParameterSetDescription &desc)
void addDefault(ParameterSetDescription const &psetDescription)
std::vector< T > JetCollection
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
def move
Definition: eostools.py:510
const double infinity
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:24
ParameterDescriptionBase * add(U const &iLabel, T const &value)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
edm::EDGetTokenT< double > m_rho_token
static const double MIN_JET_ENERGY
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
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
std::unique_ptr< JME::JetResolutionScaleFactor > m_scale_factor_from_file
tuple cout
Definition: gather_cfg.py:145
edm::EDGetTokenT< reco::GenJetCollection > m_genJetsToken
virtual double eta() const final
momentum pseudorapidity
volatile std::atomic< bool > shutdown_flag false
long double T
static std::string const source
Definition: EdmProvDump.cc:43
virtual double pt() const final
transverse momentum