CMS 3D CMS Logo

PileupJetIdProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PileupJetIdProducer
4 // Class: PileupJetIdProducer
5 //
13 //
14 // Original Author: Pasquale Musella,40 2-A12,+41227671706,
15 // Created: Wed Apr 18 15:48:47 CEST 2012
16 //
17 //
18 
20 
22  : runMvas_(iConfig.getParameter<bool>("runMvas")),
23  produceJetIds_(iConfig.getParameter<bool>("produceJetIds")),
24  inputIsCorrected_(iConfig.getParameter<bool>("inputIsCorrected")),
25  applyJec_(iConfig.getParameter<bool>("applyJec")),
26  jec_(iConfig.getParameter<std::string>("jec")),
27  residualsFromTxt_(iConfig.getParameter<bool>("residualsFromTxt")),
28  usePuppi_(iConfig.getParameter<bool>("usePuppi")) {
29  if (residualsFromTxt_) {
30  residualsTxt_ = iConfig.getParameter<edm::FileInPath>("residualsTxt");
31  }
32 
33  std::vector<edm::ParameterSet> algos = iConfig.getParameter<std::vector<edm::ParameterSet>>("algos");
34  for (auto const& algoPset : algos) {
35  vAlgoGBRForestsAndConstants_.emplace_back(algoPset, runMvas_);
36  }
37 
38  if (!runMvas_) {
39  assert(algos.size() == 1);
40  }
41 }
42 
43 // ------------------------------------------------------------------------------------------
45  if (globalCache->produceJetIds()) {
46  produces<edm::ValueMap<StoredPileupJetIdentifier>>("");
47  }
48  for (auto const& algoGBRForestsAndConstants : globalCache->vAlgoGBRForestsAndConstants()) {
49  std::string const& label = algoGBRForestsAndConstants.label();
50  algos_.emplace_back(label, std::make_unique<PileupJetIdAlgo>(&algoGBRForestsAndConstants));
51  if (globalCache->runMvas()) {
52  produces<edm::ValueMap<float>>(label + "Discriminant");
53  produces<edm::ValueMap<int>>(label + "Id");
54  }
55  }
56 
57  input_jet_token_ = consumes<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("jets"));
58  input_vertex_token_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexes"));
59  input_vm_pujetid_token_ =
60  consumes<edm::ValueMap<StoredPileupJetIdentifier>>(iConfig.getParameter<edm::InputTag>("jetids"));
61  input_rho_token_ = consumes<double>(iConfig.getParameter<edm::InputTag>("rho"));
62 }
63 
64 // ------------------------------------------------------------------------------------------
66 
67 // ------------------------------------------------------------------------------------------
69  GBRForestsAndConstants const* gc = globalCache();
70 
71  using namespace edm;
72  using namespace std;
73  using namespace reco;
74 
75  // Input jets
76  Handle<View<Jet>> jetHandle;
77  iEvent.getByToken(input_jet_token_, jetHandle);
78  const View<Jet>& jets = *jetHandle;
79 
80  // input variables
82  if (!gc->produceJetIds()) {
83  iEvent.getByToken(input_vm_pujetid_token_, vmap);
84  }
85  // rho
87  double rho = 0.;
88 
89  // products
90  vector<StoredPileupJetIdentifier> ids;
91  map<string, vector<float>> mvas;
92  map<string, vector<int>> idflags;
93 
94  const VertexCollection* vertexes = nullptr;
95  VertexCollection::const_iterator vtx;
96  if (gc->produceJetIds()) {
97  // vertexes
98  Handle<VertexCollection> vertexHandle;
99  iEvent.getByToken(input_vertex_token_, vertexHandle);
100 
101  vertexes = vertexHandle.product();
102 
103  // require basic quality cuts on the vertexes
104  vtx = vertexes->begin();
105  while (vtx != vertexes->end() && (vtx->isFake() || vtx->ndof() < 4)) {
106  ++vtx;
107  }
108  if (vtx == vertexes->end()) {
109  vtx = vertexes->begin();
110  }
111  }
112 
113  // Loop over input jets
114  bool ispat = true;
115  for (unsigned int i = 0; i < jets.size(); ++i) {
116  // Pick the first algo to compute the input variables
117  auto algoi = algos_.begin();
118  PileupJetIdAlgo* ialgo = algoi->second.get();
119 
120  const Jet& jet = jets.at(i);
121  const pat::Jet* patjet = nullptr;
122  if (ispat) {
123  patjet = dynamic_cast<const pat::Jet*>(&jet);
124  ispat = patjet != nullptr;
125  }
126 
127  // Get jet energy correction
128  float jec = 0.;
129  if (gc->applyJec()) {
130  // If haven't done it get rho from the event
131  if (rho == 0.) {
132  iEvent.getByToken(input_rho_token_, rhoH);
133  rho = *rhoH;
134  }
135  // jet corrector
136  if (not jecCor_) {
137  initJetEnergyCorrector(iSetup, iEvent.isRealData());
138  }
139  if (ispat) {
140  jecCor_->setJetPt(patjet->correctedJet(0).pt());
141  } else {
142  jecCor_->setJetPt(jet.pt());
143  }
144  jecCor_->setJetEta(jet.eta());
145  jecCor_->setJetA(jet.jetArea());
146  jecCor_->setRho(rho);
147  jec = jecCor_->getCorrection();
148  }
149  // If it was requested AND the input is an uncorrected jet apply the JEC
150  bool applyJec = gc->applyJec() && (ispat || !gc->inputIsCorrected());
151  std::unique_ptr<reco::Jet> corrJet;
152 
153  if (applyJec) {
154  float scale = jec;
155  if (ispat) {
156  corrJet.reset(new pat::Jet(patjet->correctedJet(0)));
157  } else {
158  corrJet.reset(dynamic_cast<reco::Jet*>(jet.clone()));
159  }
160  corrJet->scaleEnergy(scale);
161  }
162  const reco::Jet* theJet = (applyJec ? corrJet.get() : &jet);
163 
164  PileupJetIdentifier puIdentifier;
165  if (gc->produceJetIds()) {
166  // Compute the input variables
167  puIdentifier = ialgo->computeIdVariables(theJet, jec, &(*vtx), *vertexes, rho, gc->usePuppi());
168  ids.push_back(puIdentifier);
169  } else {
170  // Or read it from the value map
171  puIdentifier = (*vmap)[jets.refAt(i)];
172  puIdentifier.jetPt(theJet->pt()); // make sure JEC is applied when computing the MVA
173  puIdentifier.jetEta(theJet->eta());
174  puIdentifier.jetPhi(theJet->phi());
175  ialgo->set(puIdentifier);
176  puIdentifier = ialgo->computeMva();
177  }
178 
179  if (gc->runMvas()) {
180  // Compute the MVA and WP
181  mvas[algoi->first].push_back(puIdentifier.mva());
182  idflags[algoi->first].push_back(puIdentifier.idFlag());
183  for (++algoi; algoi != algos_.end(); ++algoi) {
184  ialgo = algoi->second.get();
185  ialgo->set(puIdentifier);
186  PileupJetIdentifier id = ialgo->computeMva();
187  mvas[algoi->first].push_back(id.mva());
188  idflags[algoi->first].push_back(id.idFlag());
189  }
190  }
191  }
192 
193  // Produce the output value maps
194  if (gc->runMvas()) {
195  for (const auto& ialgo : algos_) {
196  // MVA
197  vector<float>& mva = mvas[ialgo.first];
198  auto mvaout = std::make_unique<ValueMap<float>>();
199  ValueMap<float>::Filler mvafiller(*mvaout);
200  mvafiller.insert(jetHandle, mva.begin(), mva.end());
201  mvafiller.fill();
202  iEvent.put(std::move(mvaout), ialgo.first + "Discriminant");
203 
204  // WP
205  vector<int>& idflag = idflags[ialgo.first];
206  auto idflagout = std::make_unique<ValueMap<int>>();
207  ValueMap<int>::Filler idflagfiller(*idflagout);
208  idflagfiller.insert(jetHandle, idflag.begin(), idflag.end());
209  idflagfiller.fill();
210  iEvent.put(std::move(idflagout), ialgo.first + "Id");
211  }
212  }
213  // input variables
214  if (gc->produceJetIds()) {
215  assert(jetHandle->size() == ids.size());
216  auto idsout = std::make_unique<ValueMap<StoredPileupJetIdentifier>>();
218  idsfiller.insert(jetHandle, ids.begin(), ids.end());
219  idsfiller.fill();
220  iEvent.put(std::move(idsout));
221  }
222 }
223 
224 // ------------------------------------------------------------------------------------------
226  //The following says we do not know what parameters are allowed so do no validation
227  // Please change this to state exactly what you do use, even if it is no parameters
229  desc.setUnknown();
230  descriptions.addDefault(desc);
231 }
232 
233 // ------------------------------------------------------------------------------------------
235  GBRForestsAndConstants const* gc = globalCache();
236 
237  //jet energy correction levels to apply on raw jet
238  std::vector<std::string> jecLevels;
239  jecLevels.push_back("L1FastJet");
240  jecLevels.push_back("L2Relative");
241  jecLevels.push_back("L3Absolute");
242  if (isData && !gc->residualsFromTxt())
243  jecLevels.push_back("L2L3Residual");
244 
245  //check the corrector parameters needed according to the correction levels
247  iSetup.get<JetCorrectionsRecord>().get(gc->jec(), parameters);
248  for (std::vector<std::string>::const_iterator ll = jecLevels.begin(); ll != jecLevels.end(); ++ll) {
249  const JetCorrectorParameters& ip = (*parameters)[*ll];
250  jetCorPars_.push_back(ip);
251  }
252  if (isData && gc->residualsFromTxt()) {
253  jetCorPars_.push_back(JetCorrectorParameters(gc->residualsTxt().fullPath()));
254  }
255 
256  //instantiate the jet corrector
257  jecCor_ = std::make_unique<FactorizedJetCorrector>(jetCorPars_);
258 }
259 //define this as a plug-in
std::vector< PileupJetIdAlgo::AlgoGBRForestsAndConstants > const & vAlgoGBRForestsAndConstants() const
T getParameter(std::string const &) const
void set(const PileupJetIdentifier &)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const int & idFlag() const
const float & mva() const
Base class for all types of Jets.
Definition: Jet.h:20
size_type size() const
double pt() const final
transverse momentum
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool isRealData() const
Definition: EventBase.h:62
const float & jetPt() const
PileupJetIdProducer(const edm::ParameterSet &, GBRForestsAndConstants const *)
RefToBase< value_type > refAt(size_type i) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
char const * label
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho, bool usePuppi)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
Definition: Jet.py:1
std::string const & jec() const
edm::FileInPath const & residualsTxt() const
void initJetEnergyCorrector(const edm::EventSetup &iSetup, bool isData)
void produce(edm::Event &, const edm::EventSetup &) override
std::vector< PileupJetIdAlgo::AlgoGBRForestsAndConstants > vAlgoGBRForestsAndConstants_
const float & jetPhi() const
edm::FileInPath residualsTxt_
T const * product() const
Definition: Handle.h:69
PileupJetIdentifier computeMva()
Analysis-level calorimeter jet class.
Definition: Jet.h:77
CompositePtrCandidate * clone() const override
returns a clone of the candidate
fixed size matrix
HLT enums.
GBRForestsAndConstants(edm::ParameterSet const &)
virtual float jetArea() const
get jet area
Definition: Jet.h:103
T get() const
Definition: EventSetup.h:73
const_reference at(size_type pos) const
std::string fullPath() const
Definition: FileInPath.cc:163
Jet correctedJet(const std::string &level, const std::string &flavor="none", const std::string &set="") const
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:511
const float & jetEta() const