CMS 3D CMS Logo

L1PrefiringWeightProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ProdTutorial/L1PrefiringWeightProducer
4 // Class: L1PrefiringWeightProducer
5 //
13 //
14 // Original Author: localusers user
15 // Created: Thu, 08 Nov 2018 16:16:00 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
24 
27 
30 
34 
35 #include "TFile.h"
36 #include "TF1.h"
37 #include "TH2.h"
38 
39 #include <iostream>
41 
43 public:
45  ~L1PrefiringWeightProducer() override;
46 
47  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
48 
49 private:
50  void produce(edm::Event&, const edm::EventSetup&) override;
51 
52  double getPrefiringRateEcal(double eta, double pt, TH2F* h_prefmap, fluctuations fluctuation) const;
53  double getPrefiringRateMuon(double eta, double phi, double pt, fluctuations fluctuation) const;
54 
58 
62 
66 
74 
75  std::unique_ptr<TFile> file_prefiringmaps_;
76  std::unique_ptr<TFile> file_prefiringparams_;
77 
90 
95  const bool useEMpt_;
98  const double jetMaxMuonFraction_;
101 };
102 
104  : photons_token_(consumes<std::vector<pat::Photon> >(iConfig.getParameter<edm::InputTag>("ThePhotons"))),
105  jets_token_(consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("TheJets"))),
106  muon_token_(consumes<std::vector<pat::Muon> >(iConfig.getParameter<edm::InputTag>("TheMuons"))),
107  nonPrefiringProbToken_(produces<float>("nonPrefiringProb")),
108  nonPrefiringProbUpToken_(produces<float>("nonPrefiringProbUp")),
109  nonPrefiringProbDownToken_(produces<float>("nonPrefiringProbDown")),
110  nonPrefiringProbECALToken_(produces<float>("nonPrefiringProbECAL")),
111  nonPrefiringProbECALUpToken_(produces<float>("nonPrefiringProbECALUp")),
112  nonPrefiringProbECALDownToken_(produces<float>("nonPrefiringProbECALDown")),
113  nonPrefiringProbMuonToken_(produces<float>("nonPrefiringProbMuon")),
114  nonPrefiringProbMuonUpToken_(produces<float>("nonPrefiringProbMuonUp")),
115  nonPrefiringProbMuonDownToken_(produces<float>("nonPrefiringProbMuonDown")),
116  nonPrefiringProbMuonUpSystToken_(produces<float>("nonPrefiringProbMuonSystUp")),
117  nonPrefiringProbMuonDownSystToken_(produces<float>("nonPrefiringProbMuonSystDown")),
118  nonPrefiringProbMuonUpStatToken_(produces<float>("nonPrefiringProbMuonStatUp")),
119  nonPrefiringProbMuonDownStatToken_(produces<float>("nonPrefiringProbMuonStatDown")),
120  dataeraEcal_(iConfig.getParameter<std::string>("DataEraECAL")),
121  dataeraMuon_(iConfig.getParameter<std::string>("DataEraMuon")),
122  useEMpt_(iConfig.getParameter<bool>("UseJetEMPt")),
123  prefiringRateSystUncEcal_(iConfig.getParameter<double>("PrefiringRateSystematicUnctyECAL")),
124  prefiringRateSystUncMuon_(iConfig.getParameter<double>("PrefiringRateSystematicUnctyMuon")),
125  jetMaxMuonFraction_(iConfig.getParameter<double>("JetMaxMuonFraction")) {
126  missingInputEcal_ = false;
127  missingInputMuon_ = false;
128 
129  std::string fname = iConfig.getParameter<std::string>("L1Maps");
130  edm::FileInPath mapsfilepath("PhysicsTools/PatUtils/data/" + fname);
131  file_prefiringmaps_ = std::make_unique<TFile>(mapsfilepath.fullPath().c_str(), "read");
132  if (file_prefiringmaps_ == nullptr) {
133  missingInputEcal_ = true;
134  edm::LogError("L1PrefireWeightProducer")
135  << "File with maps not found. All prefiring weights set to 0. " << std::endl;
136  }
137  TString mapphotonfullname = "L1prefiring_photonptvseta_" + dataeraEcal_;
138  if (!file_prefiringmaps_->Get(mapphotonfullname)) {
139  missingInputEcal_ = true;
140  edm::LogError("L1PrefireWeightProducer")
141  << "Photon map not found. All photons prefiring weights set to 0. " << std::endl;
142  }
143  h_prefmap_photon_ = file_prefiringmaps_->Get<TH2F>(mapphotonfullname);
144  TString mapjetfullname =
145  (useEMpt_) ? "L1prefiring_jetemptvseta_" + dataeraEcal_ : "L1prefiring_jetptvseta_" + dataeraEcal_;
146  if (!file_prefiringmaps_->Get(mapjetfullname)) {
147  missingInputEcal_ = true;
148  edm::LogError("L1PrefireWeightProducer") << "Jet map not found. All jets prefiring weights set to 0. " << std::endl;
149  }
150  h_prefmap_jet_ = file_prefiringmaps_->Get<TH2F>(mapjetfullname);
151  file_prefiringmaps_->Close();
152 
153  std::string fnameMuon = iConfig.getParameter<std::string>("L1MuonParametrizations");
154  edm::FileInPath paramsfilepath("PhysicsTools/PatUtils/data/" + fnameMuon);
155  file_prefiringparams_ = std::make_unique<TFile>(paramsfilepath.fullPath().c_str(), "read");
156  if (file_prefiringparams_ == nullptr) {
157  missingInputMuon_ = true;
158  edm::LogError("L1PrefireWeightProducer")
159  << "File with muon parametrizations not found. All prefiring weights set to 0." << std::endl;
160  }
161  TString paramName = "L1prefiring_muonparam_0.0To0.2_" + dataeraMuon_;
162  parametrization0p0To0p2_ = file_prefiringparams_->Get<TF1>(paramName);
163  paramName = "L1prefiring_muonparam_0.2To0.3_" + dataeraMuon_;
164  parametrization0p2To0p3_ = file_prefiringparams_->Get<TF1>(paramName);
165  paramName = "L1prefiring_muonparam_0.3To0.55_" + dataeraMuon_;
166  parametrization0p3To0p55_ = file_prefiringparams_->Get<TF1>(paramName);
167  paramName = "L1prefiring_muonparam_0.55To0.83_" + dataeraMuon_;
168  parametrization0p55To0p83_ = file_prefiringparams_->Get<TF1>(paramName);
169  paramName = "L1prefiring_muonparam_0.83To1.24_" + dataeraMuon_;
170  parametrization0p83To1p24_ = file_prefiringparams_->Get<TF1>(paramName);
171  paramName = "L1prefiring_muonparam_1.24To1.4_" + dataeraMuon_;
172  parametrization1p24To1p4_ = file_prefiringparams_->Get<TF1>(paramName);
173  paramName = "L1prefiring_muonparam_1.4To1.6_" + dataeraMuon_;
174  parametrization1p4To1p6_ = file_prefiringparams_->Get<TF1>(paramName);
175  paramName = "L1prefiring_muonparam_1.6To1.8_" + dataeraMuon_;
176  parametrization1p6To1p8_ = file_prefiringparams_->Get<TF1>(paramName);
177  paramName = "L1prefiring_muonparam_1.8To2.1_" + dataeraMuon_;
178  parametrization1p8To2p1_ = file_prefiringparams_->Get<TF1>(paramName);
179  paramName = "L1prefiring_muonparam_2.1To2.25_" + dataeraMuon_;
180  parametrization2p1To2p25_ = file_prefiringparams_->Get<TF1>(paramName);
181  paramName = "L1prefiring_muonparam_2.25To2.4_" + dataeraMuon_;
182  parametrization2p25To2p4_ = file_prefiringparams_->Get<TF1>(paramName);
183 
184  if (parametrization0p0To0p2_ == nullptr || parametrization0p2To0p3_ == nullptr ||
185  parametrization0p3To0p55_ == nullptr || parametrization0p55To0p83_ == nullptr ||
186  parametrization0p83To1p24_ == nullptr || parametrization1p24To1p4_ == nullptr ||
187  parametrization1p4To1p6_ == nullptr || parametrization1p6To1p8_ == nullptr ||
188  parametrization1p8To2p1_ == nullptr || parametrization2p1To2p25_ == nullptr ||
189  parametrization2p25To2p4_ == nullptr) {
190  missingInputMuon_ = true;
191  edm::LogError("L1PrefireWeightProducer")
192  << "Muon parametrization not found for at least one bin. All prefiring weights set to 0." << std::endl;
193  }
194 
195  paramName = "L1prefiring_muonparam_HotSpot_" + dataeraMuon_;
196  parametrizationHotSpot_ = file_prefiringparams_->Get<TF1>(paramName);
197  file_prefiringparams_->Close();
198  if ((dataeraMuon_.find("2016") != std::string::npos) && parametrizationHotSpot_ == nullptr) {
199  missingInputMuon_ = true;
200  edm::LogError("L1PrefireWeightProducer")
201  << "Year is 2016 and no Muon parametrization is found for hot spot. All prefiring weights set to 0."
202  << std::endl;
203  }
204 }
205 
207 
209  using namespace edm;
210 
211  //Photons
212  const std::vector<pat::Photon>& thePhotons = iEvent.get(photons_token_);
213 
214  //Jets
215  const std::vector<pat::Jet>& theJets = iEvent.get(jets_token_);
216 
217  //Muons
218  const std::vector<pat::Muon>& theMuons = iEvent.get(muon_token_);
219 
220  //Probability for the event NOT to prefire, computed with the prefiring maps per object.
221  //Up and down values correspond to the resulting value when shifting up/down all prefiring rates in prefiring maps.
222  double nonPrefiringProba[3] = {1., 1., 1.}; //0: central, 1: up, 2: down
223  double nonPrefiringProbaECAL[3] = {1., 1., 1.}; //0: central, 1: up, 2: down
224  double nonPrefiringProbaMuon[7] = {
225  1., 1., 1., 1., 1., 1., 1.}; //0: central, 1: up, 2: down, 3: up stat, 4: down stat, 5: up syst, 6: down syst
226 
227  for (const auto fluct : {fluctuations::central, fluctuations::up, fluctuations::down}) {
228  if (!missingInputEcal_) {
229  for (const auto& photon : thePhotons) {
230  double pt_gam = photon.pt();
231  double eta_gam = photon.eta();
232  if (pt_gam < 20.)
233  continue;
234  if (fabs(eta_gam) < 2.)
235  continue;
236  if (fabs(eta_gam) > 3.)
237  continue;
238  double prefiringprob_gam = getPrefiringRateEcal(eta_gam, pt_gam, h_prefmap_photon_, fluct);
239  nonPrefiringProbaECAL[fluct] *= (1. - prefiringprob_gam);
240  }
241 
242  //Now applying the prefiring maps to jets in the affected regions.
243  for (const auto& jet : theJets) {
244  double pt_jet = jet.pt();
245  double eta_jet = jet.eta();
246  double phi_jet = jet.phi();
247  if (pt_jet < 20.)
248  continue;
249  if (fabs(eta_jet) < 2.)
250  continue;
251  if (fabs(eta_jet) > 3.)
252  continue;
253  if (jetMaxMuonFraction_ > 0 && jet.muonEnergyFraction() > jetMaxMuonFraction_)
254  continue;
255  //Loop over photons to remove overlap
256  double nonprefiringprobfromoverlappingphotons = 1.;
257  bool foundOverlappingPhotons = false;
258  for (const auto& photon : thePhotons) {
259  double pt_gam = photon.pt();
260  double eta_gam = photon.eta();
261  double phi_gam = photon.phi();
262  if (pt_gam < 20.)
263  continue;
264  if (fabs(eta_gam) < 2.)
265  continue;
266  if (fabs(eta_gam) > 3.)
267  continue;
268  double dR2 = reco::deltaR2(eta_jet, phi_jet, eta_gam, phi_gam);
269  if (dR2 > 0.16)
270  continue;
271  double prefiringprob_gam = getPrefiringRateEcal(eta_gam, pt_gam, h_prefmap_photon_, fluct);
272  nonprefiringprobfromoverlappingphotons *= (1. - prefiringprob_gam);
273  foundOverlappingPhotons = true;
274  }
275  //useEMpt =true if one wants to use maps parametrized vs Jet EM pt instead of pt.
276  if (useEMpt_)
277  pt_jet *= (jet.neutralEmEnergyFraction() + jet.chargedEmEnergyFraction());
278  double nonprefiringprobfromoverlappingjet = 1. - getPrefiringRateEcal(eta_jet, pt_jet, h_prefmap_jet_, fluct);
279 
280  if (!foundOverlappingPhotons) {
281  nonPrefiringProbaECAL[fluct] *= nonprefiringprobfromoverlappingjet;
282  }
283  //If overlapping photons have a non prefiring rate larger than the jet, then replace these weights by the jet one
284  else if (nonprefiringprobfromoverlappingphotons > nonprefiringprobfromoverlappingjet) {
285  if (nonprefiringprobfromoverlappingphotons > 0.) {
286  nonPrefiringProbaECAL[fluct] *= nonprefiringprobfromoverlappingjet / nonprefiringprobfromoverlappingphotons;
287  } else {
288  nonPrefiringProbaECAL[fluct] = 0.;
289  }
290  }
291  //Last case: if overlapping photons have a non prefiring rate smaller than the jet, don't consider the jet in the event weight, and do nothing.
292  }
293  }
294  //Now calculate prefiring weights for muons
295  if (!missingInputMuon_) {
296  for (const auto& muon : theMuons) {
297  double pt = muon.pt();
298  double phi = muon.phi();
299  double eta = muon.eta();
300  // Remove crappy tracker muons which would not have prefired the L1 trigger
301  if (pt < 5 || !muon.isLooseMuon())
302  continue;
303  double prefiringprob_mu = getPrefiringRateMuon(eta, phi, pt, fluct);
304  nonPrefiringProbaMuon[fluct] *= (1. - prefiringprob_mu);
305  }
306  }
307  }
308  // Calculate combined weight as product of the weight for individual objects
309  for (const auto fluct : {fluctuations::central, fluctuations::up, fluctuations::down}) {
310  nonPrefiringProba[fluct] = nonPrefiringProbaECAL[fluct] * nonPrefiringProbaMuon[fluct];
311  }
312  // Calculate statistical and systematic uncertainty separately in the muon case
313  for (const auto fluct :
315  if (!missingInputMuon_) {
316  for (const auto& muon : theMuons) {
317  double pt = muon.pt();
318  double phi = muon.phi();
319  double eta = muon.eta();
320  // Remove crappy tracker muons which would not have prefired the L1 trigger
321  if (pt < 5 || !muon.isLooseMuon())
322  continue;
323  double prefiringprob_mu = getPrefiringRateMuon(eta, phi, pt, fluct);
324  nonPrefiringProbaMuon[fluct] *= (1. - prefiringprob_mu);
325  }
326  }
327  }
328  //Move global prefire weights, as well as those for muons, photons, and jets, to the event
329  iEvent.emplace(nonPrefiringProbToken_, nonPrefiringProba[0]);
330  iEvent.emplace(nonPrefiringProbUpToken_, nonPrefiringProba[1]);
331  iEvent.emplace(nonPrefiringProbDownToken_, nonPrefiringProba[2]);
332 
333  iEvent.emplace(nonPrefiringProbECALToken_, nonPrefiringProbaECAL[0]);
334  iEvent.emplace(nonPrefiringProbECALUpToken_, nonPrefiringProbaECAL[1]);
335  iEvent.emplace(nonPrefiringProbECALDownToken_, nonPrefiringProbaECAL[2]);
336 
337  iEvent.emplace(nonPrefiringProbMuonToken_, nonPrefiringProbaMuon[0]);
338  iEvent.emplace(nonPrefiringProbMuonUpToken_, nonPrefiringProbaMuon[1]);
339  iEvent.emplace(nonPrefiringProbMuonDownToken_, nonPrefiringProbaMuon[2]);
340  iEvent.emplace(nonPrefiringProbMuonUpStatToken_, nonPrefiringProbaMuon[3]);
341  iEvent.emplace(nonPrefiringProbMuonDownStatToken_, nonPrefiringProbaMuon[4]);
342  iEvent.emplace(nonPrefiringProbMuonUpSystToken_, nonPrefiringProbaMuon[5]);
343  iEvent.emplace(nonPrefiringProbMuonDownSystToken_, nonPrefiringProbaMuon[6]);
344 }
345 
347  double pt,
348  TH2F* h_prefmap,
349  fluctuations fluctuation) const {
350  //Check pt is not above map overflow
351  int nbinsy = h_prefmap->GetNbinsY();
352  double maxy = h_prefmap->GetYaxis()->GetBinLowEdge(nbinsy + 1);
353  if (pt >= maxy)
354  pt = maxy - 0.01;
355  int thebin = h_prefmap->FindBin(eta, pt);
356 
357  double prefrate = h_prefmap->GetBinContent(thebin);
358 
359  double statuncty = h_prefmap->GetBinError(thebin);
360  double systuncty = prefiringRateSystUncEcal_ * prefrate;
361 
362  if (fluctuation == up)
363  prefrate = std::min(1., prefrate + sqrt(pow(statuncty, 2) + pow(systuncty, 2)));
364  else if (fluctuation == down)
365  prefrate = std::max(0., prefrate - sqrt(pow(statuncty, 2) + pow(systuncty, 2)));
366  if (prefrate > 1.) {
367  edm::LogWarning("L1PrefireWeightProducer") << "Found a prefiring probability > 1. Setting to 1." << std::endl;
368  return 1.;
369  }
370  return prefrate;
371 }
372 
374  double phi,
375  double pt,
376  fluctuations fluctuation) const {
377  double prefrate;
378  double statuncty;
379  if ((dataeraMuon_.find("2016") != std::string::npos) && (eta > 1.24 && eta < 1.6) &&
380  (phi > 2.44346 && phi < 2.79253)) {
381  prefrate = parametrizationHotSpot_->Eval(pt);
382  statuncty = parametrizationHotSpot_->GetParError(2);
383  } else if (std::abs(eta) < 0.2) {
384  prefrate = parametrization0p0To0p2_->Eval(pt);
385  statuncty = parametrization0p0To0p2_->GetParError(2);
386  } else if (std::abs(eta) < 0.3) {
387  prefrate = parametrization0p2To0p3_->Eval(pt);
388  statuncty = parametrization0p2To0p3_->GetParError(2);
389  } else if (std::abs(eta) < 0.55) {
390  prefrate = parametrization0p3To0p55_->Eval(pt);
391  statuncty = parametrization0p3To0p55_->GetParError(2);
392  } else if (std::abs(eta) < 0.83) {
393  prefrate = parametrization0p55To0p83_->Eval(pt);
394  statuncty = parametrization0p55To0p83_->GetParError(2);
395  } else if (std::abs(eta) < 1.24) {
396  prefrate = parametrization0p83To1p24_->Eval(pt);
397  statuncty = parametrization0p83To1p24_->GetParError(2);
398  } else if (std::abs(eta) < 1.4) {
399  prefrate = parametrization1p24To1p4_->Eval(pt);
400  statuncty = parametrization1p24To1p4_->GetParError(2);
401  } else if (std::abs(eta) < 1.6) {
402  prefrate = parametrization1p4To1p6_->Eval(pt);
403  statuncty = parametrization1p4To1p6_->GetParError(2);
404  } else if (std::abs(eta) < 1.8) {
405  prefrate = parametrization1p6To1p8_->Eval(pt);
406  statuncty = parametrization1p6To1p8_->GetParError(2);
407  } else if (std::abs(eta) < 2.1) {
408  prefrate = parametrization1p8To2p1_->Eval(pt);
409  statuncty = parametrization1p8To2p1_->GetParError(2);
410  } else if (std::abs(eta) < 2.25) {
411  prefrate = parametrization2p1To2p25_->Eval(pt);
412  statuncty = parametrization2p1To2p25_->GetParError(2);
413  } else if (std::abs(eta) < 2.4) {
414  prefrate = parametrization2p25To2p4_->Eval(pt);
415  statuncty = parametrization2p25To2p4_->GetParError(2);
416  } else {
417  LogDebug("L1PrefireWeightProducer") << "Muon outside of |eta| <= 2.4. Prefiring weight set to 0." << std::endl;
418  return 0.;
419  }
420  double systuncty = prefiringRateSystUncMuon_ * prefrate;
421 
422  if (fluctuation == up)
423  prefrate = std::min(1., prefrate + sqrt(pow(statuncty, 2) + pow(systuncty, 2)));
424  else if (fluctuation == down)
425  prefrate = std::max(0., prefrate - sqrt(pow(statuncty, 2) + pow(systuncty, 2)));
426  else if (fluctuation == upSyst)
427  prefrate = std::min(1., prefrate + systuncty);
428  else if (fluctuation == downSyst)
429  prefrate = std::max(0., prefrate - systuncty);
430  else if (fluctuation == upStat)
431  prefrate = std::min(1., prefrate + statuncty);
432  else if (fluctuation == downStat)
433  prefrate = std::max(0., prefrate - statuncty);
434 
435  if (prefrate > 1.) {
436  edm::LogWarning("L1PrefireWeightProducer") << "Found a prefiring probability > 1. Setting to 1." << std::endl;
437  return 1.;
438  }
439  return prefrate;
440 }
441 
442 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
445  desc.add<edm::InputTag>("TheMuons", edm::InputTag("slimmedMuons"));
446  desc.add<edm::InputTag>("ThePhotons", edm::InputTag("slimmedPhotons"));
447  desc.add<edm::InputTag>("TheJets", edm::InputTag("slimmedJets"));
448  desc.add<std::string>("L1Maps", "L1PrefiringMaps.root");
449  desc.add<std::string>("L1MuonParametrizations", "L1MuonPrefiringParametriations.root");
450  desc.add<std::string>("DataEraECAL", "2017BtoF");
451  desc.add<std::string>("DataEraMuon", "2016");
452  desc.add<bool>("UseJetEMPt", false);
453  desc.add<double>("PrefiringRateSystematicUnctyECAL", 0.2);
454  desc.add<double>("PrefiringRateSystematicUnctyMuon", 0.2);
455  desc.add<double>("JetMaxMuonFraction", 0.5);
456  descriptions.add("l1PrefiringWeightProducer", desc);
457 }
458 
459 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Definition: Photon.py:1
const edm::EDPutTokenT< float > nonPrefiringProbDownToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDPutTokenT< float > nonPrefiringProbECALUpToken_
L1PrefiringWeightProducer(const edm::ParameterSet &)
const edm::EDPutTokenT< float > nonPrefiringProbMuonDownSystToken_
const edm::EDPutTokenT< float > nonPrefiringProbMuonToken_
constexpr int pow(int x)
Definition: conifer.h:24
Log< level::Error, false > LogError
std::unique_ptr< TFile > file_prefiringparams_
const edm::EDGetTokenT< std::vector< pat::Photon > > photons_token_
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
double getPrefiringRateMuon(double eta, double phi, double pt, fluctuations fluctuation) const
Definition: Muon.py:1
const edm::EDPutTokenT< float > nonPrefiringProbECALToken_
Definition: Jet.py:1
T sqrt(T t)
Definition: SSEVec.h:19
const edm::EDPutTokenT< float > nonPrefiringProbECALDownToken_
void produce(edm::Event &, const edm::EventSetup &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDPutTokenT< float > nonPrefiringProbMuonUpSystToken_
const edm::EDPutTokenT< float > nonPrefiringProbMuonUpStatToken_
const edm::EDPutTokenT< float > nonPrefiringProbMuonUpToken_
const edm::EDPutTokenT< float > nonPrefiringProbToken_
const edm::EDPutTokenT< float > nonPrefiringProbUpToken_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
double getPrefiringRateEcal(double eta, double pt, TH2F *h_prefmap, fluctuations fluctuation) const
const edm::EDPutTokenT< float > nonPrefiringProbMuonDownToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
string fname
main script
const edm::EDPutTokenT< float > nonPrefiringProbMuonDownStatToken_
HLT enums.
const edm::EDGetTokenT< std::vector< pat::Muon > > muon_token_
const edm::EDGetTokenT< std::vector< pat::Jet > > jets_token_
Log< level::Warning, false > LogWarning
std::unique_ptr< TFile > file_prefiringmaps_
#define LogDebug(id)