CMS 3D CMS Logo

L1ECALPrefiringWeightProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ProdTutorial/L1ECALPrefiringWeightProducer
4 // Class: L1ECALPrefiringWeightProducer
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 
33 #include "TH2.h"
34 
35 #include <iostream>
36 enum fluctuations { central = 0, up, down };
37 
39 public:
42 
43  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
44 
45 private:
46  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
47 
48  double getPrefiringRate(double eta, double pt, TH2F* h_prefmap, fluctuations fluctuation) const;
49 
52 
56  bool useEMpt_;
60 };
61 
63  photons_token_ = consumes<std::vector<pat::Photon> >(iConfig.getParameter<edm::InputTag>("ThePhotons"));
64  jets_token_ = consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("TheJets"));
65 
66  dataera_ = iConfig.getParameter<std::string>("DataEra");
67  useEMpt_ = iConfig.getParameter<bool>("UseJetEMPt");
68  prefiringRateSystUnc_ = iConfig.getParameter<double>("PrefiringRateSystematicUncty");
69  jetMaxMuonFraction_ = iConfig.getParameter<double>("JetMaxMuonFraction");
70  skipwarnings_ = iConfig.getParameter<bool>("SkipWarnings");
71 
72  TFile* file_prefiringmaps_;
73  std::string fname = iConfig.getParameter<std::string>("L1Maps");
74  edm::FileInPath mapsfilepath("PhysicsTools/PatUtils/data/" + fname);
75  file_prefiringmaps_ = new TFile(mapsfilepath.fullPath().c_str(), "read");
76  if (file_prefiringmaps_ == nullptr && !skipwarnings_)
77  std::cout << "File with maps not found. All prefiring weights set to 0. " << std::endl;
78 
79  TString mapphotonfullname = "L1prefiring_photonptvseta_" + dataera_;
80  if (!file_prefiringmaps_->Get(mapphotonfullname) && !skipwarnings_)
81  std::cout << "Photon map not found. All photons prefiring weights set to 0. " << std::endl;
82  h_prefmap_photon = (TH2F*)file_prefiringmaps_->Get(mapphotonfullname);
83 
84  TString mapjetfullname = (useEMpt_) ? "L1prefiring_jetemptvseta_" + dataera_ : "L1prefiring_jetptvseta_" + dataera_;
85  if (!file_prefiringmaps_->Get(mapjetfullname) && !skipwarnings_)
86  std::cout << "Jet map not found. All jets prefiring weights set to 0. " << std::endl;
87  h_prefmap_jet = (TH2F*)file_prefiringmaps_->Get(mapjetfullname);
88  file_prefiringmaps_->Close();
89  delete file_prefiringmaps_;
90  produces<double>("nonPrefiringProb").setBranchAlias("nonPrefiringProb");
91  produces<double>("nonPrefiringProbUp").setBranchAlias("nonPrefiringProbUp");
92  produces<double>("nonPrefiringProbDown").setBranchAlias("nonPrefiringProbDown");
93 }
94 
96  delete h_prefmap_photon;
97  delete h_prefmap_jet;
98 }
99 
101  using namespace edm;
102 
103  //Photons
105  iEvent.getByToken(photons_token_, thePhotons);
106 
107  //Jets
109  iEvent.getByToken(jets_token_, theJets);
110 
111  //Probability for the event NOT to prefire, computed with the prefiring maps per object.
112  //Up and down values correspond to the resulting value when shifting up/down all prefiring rates in prefiring maps.
113  double nonPrefiringProba[3] = {1., 1., 1.}; //0: central, 1: up, 2: down
114 
115  for (const auto fluct : {fluctuations::central, fluctuations::up, fluctuations::down}) {
116  for (const auto& photon : *thePhotons) {
117  double pt_gam = photon.pt();
118  double eta_gam = photon.eta();
119  if (pt_gam < 20.)
120  continue;
121  if (fabs(eta_gam) < 2.)
122  continue;
123  if (fabs(eta_gam) > 3.)
124  continue;
125  double prefiringprob_gam = getPrefiringRate(eta_gam, pt_gam, h_prefmap_photon, fluct);
126  nonPrefiringProba[fluct] *= (1. - prefiringprob_gam);
127  }
128 
129  //Now applying the prefiring maps to jets in the affected regions.
130  for (const auto& jet : *theJets) {
131  double pt_jet = jet.pt();
132  double eta_jet = jet.eta();
133  double phi_jet = jet.phi();
134  if (pt_jet < 20.)
135  continue;
136  if (fabs(eta_jet) < 2.)
137  continue;
138  if (fabs(eta_jet) > 3.)
139  continue;
140  if (jetMaxMuonFraction_ > 0 && jet.muonEnergyFraction() > jetMaxMuonFraction_)
141  continue;
142  //Loop over photons to remove overlap
143  double nonprefiringprobfromoverlappingphotons = 1.;
144  for (const auto& photon : *thePhotons) {
145  double pt_gam = photon.pt();
146  double eta_gam = photon.eta();
147  double phi_gam = photon.phi();
148  if (pt_gam < 20.)
149  continue;
150  if (fabs(eta_gam) < 2.)
151  continue;
152  if (fabs(eta_gam) > 3.)
153  continue;
154  double dR = reco::deltaR(eta_jet, phi_jet, eta_gam, phi_gam);
155  if (dR > 0.4)
156  continue;
157  double prefiringprob_gam = getPrefiringRate(eta_gam, pt_gam, h_prefmap_photon, fluct);
158  nonprefiringprobfromoverlappingphotons *= (1. - prefiringprob_gam);
159  }
160  //useEMpt =true if one wants to use maps parametrized vs Jet EM pt instead of pt.
161  if (useEMpt_)
162  pt_jet *= (jet.neutralEmEnergyFraction() + jet.chargedEmEnergyFraction());
163  double nonprefiringprobfromoverlappingjet = 1. - getPrefiringRate(eta_jet, pt_jet, h_prefmap_jet, fluct);
164 
165  if (nonprefiringprobfromoverlappingphotons == 1.) {
166  nonPrefiringProba[fluct] *= nonprefiringprobfromoverlappingjet;
167  }
168  //If overlapping photons have a non prefiring rate larger than the jet, then replace these weights by the jet one
169  else if (nonprefiringprobfromoverlappingphotons > nonprefiringprobfromoverlappingjet) {
170  if (nonprefiringprobfromoverlappingphotons != 0.) {
171  nonPrefiringProba[fluct] *= nonprefiringprobfromoverlappingjet / nonprefiringprobfromoverlappingphotons;
172  } else {
173  nonPrefiringProba[fluct] = 0.;
174  }
175  }
176  //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.
177  }
178  }
179 
180  auto nonPrefiringProb = std::make_unique<double>(nonPrefiringProba[0]);
181  auto nonPrefiringProbUp = std::make_unique<double>(nonPrefiringProba[1]);
182  auto nonPrefiringProbDown = std::make_unique<double>(nonPrefiringProba[2]);
183  iEvent.put(std::move(nonPrefiringProb), "nonPrefiringProb");
184  iEvent.put(std::move(nonPrefiringProbUp), "nonPrefiringProbUp");
185  iEvent.put(std::move(nonPrefiringProbDown), "nonPrefiringProbDown");
186 }
187 
189  double pt,
190  TH2F* h_prefmap,
191  fluctuations fluctuation) const {
192  if (h_prefmap == nullptr && !skipwarnings_)
193  std::cout << "Prefiring map not found, setting prefiring rate to 0 " << std::endl;
194  if (h_prefmap == nullptr)
195  return 0.;
196  //Check pt is not above map overflow
197  int nbinsy = h_prefmap->GetNbinsY();
198  double maxy = h_prefmap->GetYaxis()->GetBinLowEdge(nbinsy + 1);
199  if (pt >= maxy)
200  pt = maxy - 0.01;
201  int thebin = h_prefmap->FindBin(eta, pt);
202 
203  double prefrate = h_prefmap->GetBinContent(thebin);
204 
205  double statuncty = h_prefmap->GetBinError(thebin);
206  double systuncty = prefiringRateSystUnc_ * prefrate;
207 
208  if (fluctuation == up)
209  prefrate = std::min(1., prefrate + sqrt(pow(statuncty, 2) + pow(systuncty, 2)));
210  if (fluctuation == down)
211  prefrate = std::max(0., prefrate - sqrt(pow(statuncty, 2) + pow(systuncty, 2)));
212  return prefrate;
213 }
214 
215 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
218 
219  desc.add<edm::InputTag>("ThePhotons", edm::InputTag("slimmedPhotons"));
220  desc.add<edm::InputTag>("TheJets", edm::InputTag("slimmedJets"));
221  desc.add<std::string>("L1Maps", "L1PrefiringMaps.root");
222  desc.add<std::string>("DataEra", "2017BtoF");
223  desc.add<bool>("UseJetEMPt", false);
224  desc.add<double>("PrefiringRateSystematicUncty", 0.2);
225  desc.add<double>("JetMaxMuonFraction", 0.5);
226  desc.add<bool>("SkipWarnings", true);
227  descriptions.add("l1ECALPrefiringWeightProducer", desc);
228 }
229 
230 //define this as a plug-in
L1ECALPrefiringWeightProducer::photons_token_
edm::EDGetTokenT< std::vector< pat::Photon > > photons_token_
Definition: L1ECALPrefiringWeightProducer.cc:50
edm::StreamID
Definition: StreamID.h:30
L1ECALPrefiringWeightProducer::dataera_
std::string dataera_
Definition: L1ECALPrefiringWeightProducer.cc:55
muons2muons_cfi.photon
photon
Definition: muons2muons_cfi.py:28
L1ECALPrefiringWeightProducer::~L1ECALPrefiringWeightProducer
~L1ECALPrefiringWeightProducer() override
Definition: L1ECALPrefiringWeightProducer.cc:95
L1ECALPrefiringWeightProducer::skipwarnings_
bool skipwarnings_
Definition: L1ECALPrefiringWeightProducer.cc:59
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
L1ECALPrefiringWeightProducer::getPrefiringRate
double getPrefiringRate(double eta, double pt, TH2F *h_prefmap, fluctuations fluctuation) const
Definition: L1ECALPrefiringWeightProducer.cc:188
min
T min(T a, T b)
Definition: MathUtil.h:58
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
Photon.h
gather_cfg.cout
cout
Definition: gather_cfg.py:144
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89285
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
central
Definition: L1ECALPrefiringWeightProducer.cc:36
L1ECALPrefiringWeightProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L1ECALPrefiringWeightProducer.cc:216
edm::Handle
Definition: AssociativeIterator.h:50
edm::FileInPath
Definition: FileInPath.h:64
L1ECALPrefiringWeightProducer::useEMpt_
bool useEMpt_
Definition: L1ECALPrefiringWeightProducer.cc:56
MakerMacros.h
L1ECALPrefiringWeightProducer::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: L1ECALPrefiringWeightProducer.cc:100
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
PVValHelper::eta
Definition: PVValidationHelpers.h:70
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
L1ECALPrefiringWeightProducer::L1ECALPrefiringWeightProducer
L1ECALPrefiringWeightProducer(const edm::ParameterSet &)
Definition: L1ECALPrefiringWeightProducer.cc:62
L1ECALPrefiringWeightProducer::h_prefmap_jet
TH2F * h_prefmap_jet
Definition: L1ECALPrefiringWeightProducer.cc:54
edm::global::EDProducer
Definition: EDProducer.h:32
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
L1ECALPrefiringWeightProducer::h_prefmap_photon
TH2F * h_prefmap_photon
Definition: L1ECALPrefiringWeightProducer.cc:53
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
iEvent
int iEvent
Definition: GenABIO.cc:224
L1ECALPrefiringWeightProducer
Definition: L1ECALPrefiringWeightProducer.cc:38
edm::EventSetup
Definition: EventSetup.h:58
down
Definition: L1ECALPrefiringWeightProducer.cc:36
Jet.h
L1ECALPrefiringWeightProducer::jets_token_
edm::EDGetTokenT< std::vector< pat::Jet > > jets_token_
Definition: L1ECALPrefiringWeightProducer.cc:51
alignmentValidation.fname
string fname
main script
Definition: alignmentValidation.py:959
L1ECALPrefiringWeightProducer::jetMaxMuonFraction_
double jetMaxMuonFraction_
Definition: L1ECALPrefiringWeightProducer.cc:58
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
fluctuations
fluctuations
Definition: L1ECALPrefiringWeightProducer.cc:36
Frameworkfwd.h
metsig::jet
Definition: SignAlgoResolutions.h:47
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
up
Definition: L1ECALPrefiringWeightProducer.cc:36
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition: HGC3DClusterGenMatchSelector_cfi.py:7
ParameterSet.h
EDProducer.h
edm::Event
Definition: Event.h:73
StreamID.h
edm::InputTag
Definition: InputTag.h:15
L1ECALPrefiringWeightProducer::prefiringRateSystUnc_
double prefiringRateSystUnc_
Definition: L1ECALPrefiringWeightProducer.cc:57