CMS 3D CMS Logo

EGExtraInfoModifierFromDB.cc
Go to the documentation of this file.
7 
10 
14 
15 #include <vdt/vdtMath.h>
16 
17 namespace {
18  const edm::InputTag empty_tag;
19 }
20 
21 #include <unordered_map>
22 
24 public:
27  typedef std::pair<edm::InputTag, ValMapFloatToken> ValMapFloatTagTokenPair;
28  typedef std::pair<edm::InputTag, ValMapIntToken> ValMapIntTagTokenPair;
29 
30  struct electron_config {
33  std::unordered_map<std::string, ValMapFloatTagTokenPair> tag_float_token_map;
34  std::unordered_map<std::string, ValMapIntTagTokenPair> tag_int_token_map;
35 
36  std::vector<std::string> condnames_ecalonly_mean;
37  std::vector<std::string> condnames_ecalonly_sigma;
38  std::vector<std::string> condnames_ecaltrk_mean;
39  std::vector<std::string> condnames_ecaltrk_sigma;
40  };
41 
42  struct photon_config {
45  std::unordered_map<std::string, ValMapFloatTagTokenPair> tag_float_token_map;
46  std::unordered_map<std::string, ValMapIntTagTokenPair> tag_int_token_map;
47 
48  std::vector<std::string> condnames_ecalonly_mean;
49  std::vector<std::string> condnames_ecalonly_sigma;
50  };
51 
54 
55  void setEvent(const edm::Event&) override final;
56  void setEventContent(const edm::EventSetup&) override final;
57  void setConsumes(edm::ConsumesCollector&) override final;
58 
59  void modifyObject(reco::GsfElectron&) const override final;
60  void modifyObject(reco::Photon&) const override final;
61 
62  // just calls reco versions
63  void modifyObject(pat::Electron&) const override final;
64  void modifyObject(pat::Photon&) const override final;
65 
66 private:
69  std::unordered_map<unsigned,edm::Ptr<reco::GsfElectron> > eles_by_oop; // indexed by original object ptr
70  std::unordered_map<unsigned,edm::Handle<edm::ValueMap<float> > > ele_vmaps;
71  std::unordered_map<unsigned,edm::Handle<edm::ValueMap<int> > > ele_int_vmaps;
72  std::unordered_map<unsigned,edm::Ptr<reco::Photon> > phos_by_oop;
73  std::unordered_map<unsigned,edm::Handle<edm::ValueMap<float> > > pho_vmaps;
74  std::unordered_map<unsigned,edm::Handle<edm::ValueMap<int> > > pho_int_vmaps;
75 
76  float rhoValue_;
77  edm::InputTag rhoTag_;
78  edm::EDGetTokenT<double> rhoToken_;
79 
81 
82  std::vector<const GBRForestD*> ph_forestH_mean_;
83  std::vector<const GBRForestD*> ph_forestH_sigma_;
84  std::vector<const GBRForestD*> e_forestH_mean_;
85  std::vector<const GBRForestD*> e_forestH_sigma_;
86 
87  double lowEnergy_ECALonlyThr_; // 300
88  double lowEnergy_ECALTRKThr_; // 50
89  double highEnergy_ECALTRKThr_; // 200
90  double eOverP_ECALTRKThr_; // 0.025
91  double epDiffSig_ECALTRKThr_; // 15
92  double epSig_ECALTRKThr_; // 10
93 
94 
95 };
96 
99  "EGExtraInfoModifierFromDB");
100 
101 EGExtraInfoModifierFromDB::EGExtraInfoModifierFromDB(const edm::ParameterSet& conf) :
102  ModifyObjectValueBase(conf) {
103 
104  lowEnergy_ECALonlyThr_ = conf.getParameter<double>("lowEnergy_ECALonlyThr");
105  lowEnergy_ECALTRKThr_ = conf.getParameter<double>("lowEnergy_ECALTRKThr");
106  highEnergy_ECALTRKThr_ = conf.getParameter<double>("highEnergy_ECALTRKThr");
107  eOverP_ECALTRKThr_ = conf.getParameter<double>("eOverP_ECALTRKThr");
108  epDiffSig_ECALTRKThr_ = conf.getParameter<double>("epDiffSig_ECALTRKThr");
109  epSig_ECALTRKThr_ = conf.getParameter<double>("epSig_ECALTRKThr");
110 
111  rhoTag_ = conf.getParameter<edm::InputTag>("rhoCollection");
112 
113  constexpr char electronSrc[] = "electronSrc";
114  constexpr char photonSrc[] = "photonSrc";
115 
116  if(conf.exists("electron_config")) {
117  const edm::ParameterSet& electrons = conf.getParameter<edm::ParameterSet>("electron_config");
118  if( electrons.exists(electronSrc) )
119  e_conf.electron_src = electrons.getParameter<edm::InputTag>(electronSrc);
120 
121  std::vector<std::string> intValueMaps;
122  if ( electrons.existsAs<std::vector<std::string> >("intValueMaps"))
123  intValueMaps = electrons.getParameter<std::vector<std::string> >("intValueMaps");
124 
125  const std::vector<std::string> parameters = electrons.getParameterNames();
126  for( const std::string& name : parameters ) {
127  if( std::string(electronSrc) == name )
128  continue;
129  if( electrons.existsAs<edm::InputTag>(name)) {
130  for (auto vmp : intValueMaps) {
131  if (name == vmp) {
132  e_conf.tag_int_token_map[name] = ValMapIntTagTokenPair(electrons.getParameter<edm::InputTag>(name), ValMapIntToken());
133  break;
134  }
135  }
136  e_conf.tag_float_token_map[name] = ValMapFloatTagTokenPair(electrons.getParameter<edm::InputTag>(name), ValMapFloatToken());
137  }
138  }
139 
140  e_conf.condnames_ecalonly_mean = electrons.getParameter<std::vector<std::string> >("regressionKey_ecalonly");
141  e_conf.condnames_ecalonly_sigma = electrons.getParameter<std::vector<std::string> >("uncertaintyKey_ecalonly");
142  e_conf.condnames_ecaltrk_mean = electrons.getParameter<std::vector<std::string> >("regressionKey_ecaltrk");
143  e_conf.condnames_ecaltrk_sigma = electrons.getParameter<std::vector<std::string> >("uncertaintyKey_ecaltrk");
144 
145  unsigned int encor = e_conf.condnames_ecalonly_mean.size();
146  e_forestH_mean_.reserve(2*encor);
147  e_forestH_sigma_.reserve(2*encor);
148 
149  }
150 
151  if( conf.exists("photon_config") ) {
152  const edm::ParameterSet& photons = conf.getParameter<edm::ParameterSet>("photon_config");
153 
154  if( photons.exists(photonSrc) )
155  ph_conf.photon_src = photons.getParameter<edm::InputTag>(photonSrc);
156 
157  std::vector<std::string> intValueMaps;
158  if ( photons.existsAs<std::vector<std::string> >("intValueMaps"))
159  intValueMaps = photons.getParameter<std::vector<std::string> >("intValueMaps");
160 
161  const std::vector<std::string> parameters = photons.getParameterNames();
162  for( const std::string& name : parameters ) {
163  if( std::string(photonSrc) == name )
164  continue;
165  if( photons.existsAs<edm::InputTag>(name)) {
166  for (auto vmp : intValueMaps) {
167  if (name == vmp) {
168  ph_conf.tag_int_token_map[name] = ValMapIntTagTokenPair(photons.getParameter<edm::InputTag>(name), ValMapIntToken());
169  break;
170  }
171  }
172  ph_conf.tag_float_token_map[name] = ValMapFloatTagTokenPair(photons.getParameter<edm::InputTag>(name), ValMapFloatToken());
173  }
174  }
175 
176  ph_conf.condnames_ecalonly_mean = photons.getParameter<std::vector<std::string>>("regressionKey_ecalonly");
177  ph_conf.condnames_ecalonly_sigma = photons.getParameter<std::vector<std::string>>("uncertaintyKey_ecalonly");
178 
179  unsigned int ncor = ph_conf.condnames_ecalonly_mean.size();
180  ph_forestH_mean_.reserve(ncor);
181  ph_forestH_sigma_.reserve(ncor);
182 
183  }
184 
185 
186 }
187 
188 namespace {
189  template<typename T>
190  inline void get_product(const edm::Event& evt,
191  const edm::EDGetTokenT<edm::ValueMap<T> >& tok,
192  std::unordered_map<unsigned, edm::Handle<edm::ValueMap<T> > >& map) {
193  evt.getByToken(tok,map[tok.index()]);
194  }
195 }
196 
198 
200 
201  eles_by_oop.clear();
202  phos_by_oop.clear();
203  ele_vmaps.clear();
204  ele_int_vmaps.clear();
205  pho_vmaps.clear();
206  pho_int_vmaps.clear();
207 
211 
212  for( unsigned i = 0; i < eles->size(); ++i ) {
213  edm::Ptr<pat::Electron> ptr = eles->ptrAt(i);
214  eles_by_oop[ptr->originalObjectRef().key()] = ptr;
215  }
216  }
217 
218  for (std::unordered_map<std::string, ValMapFloatTagTokenPair>::iterator imap = e_conf.tag_float_token_map.begin();
219  imap != e_conf.tag_float_token_map.end();
220  imap++) {
221  get_product(evt, imap->second.second, ele_vmaps);
222  }
223 
224  for (std::unordered_map<std::string, ValMapIntTagTokenPair>::iterator imap = e_conf.tag_int_token_map.begin();
225  imap != e_conf.tag_int_token_map.end();
226  imap++) {
227  get_product(evt, imap->second.second, ele_int_vmaps);
228  }
229 
233 
234  for( unsigned i = 0; i < phos->size(); ++i ) {
235  edm::Ptr<pat::Photon> ptr = phos->ptrAt(i);
236  phos_by_oop[ptr->originalObjectRef().key()] = ptr;
237  }
238  }
239 
240 
241  for (std::unordered_map<std::string, ValMapFloatTagTokenPair>::iterator imap = ph_conf.tag_float_token_map.begin();
242  imap != ph_conf.tag_float_token_map.end();
243  imap++) {
244  get_product(evt, imap->second.second, pho_vmaps);
245  }
246 
247  for (std::unordered_map<std::string, ValMapIntTagTokenPair>::iterator imap = ph_conf.tag_int_token_map.begin();
248  imap != ph_conf.tag_int_token_map.end();
249  imap++) {
250  get_product(evt, imap->second.second, pho_int_vmaps);
251  }
252 
253  edm::Handle<double> rhoH;
254  evt.getByToken(rhoToken_, rhoH);
255  rhoValue_ = *rhoH;
256 
257 }
258 
260 
261  iSetup_ = &evs;
262 
263  edm::ESHandle<GBRForestD> forestDEH;
264 
265  const std::vector<std::string> ph_condnames_ecalonly_mean = ph_conf.condnames_ecalonly_mean;
266  const std::vector<std::string> ph_condnames_ecalonly_sigma = ph_conf.condnames_ecalonly_sigma;
267 
268  unsigned int ncor = ph_condnames_ecalonly_mean.size();
269  for (unsigned int icor=0; icor<ncor; ++icor) {
270  evs.get<GBRDWrapperRcd>().get(ph_condnames_ecalonly_mean[icor], forestDEH);
271  ph_forestH_mean_[icor] = forestDEH.product();
272  evs.get<GBRDWrapperRcd>().get(ph_condnames_ecalonly_sigma[icor], forestDEH);
273  ph_forestH_sigma_[icor] = forestDEH.product();
274  }
275 
276  const std::vector<std::string> e_condnames_ecalonly_mean = e_conf.condnames_ecalonly_mean;
277  const std::vector<std::string> e_condnames_ecalonly_sigma = e_conf.condnames_ecalonly_sigma;
278  const std::vector<std::string> e_condnames_ecaltrk_mean = e_conf.condnames_ecaltrk_mean;
279  const std::vector<std::string> e_condnames_ecaltrk_sigma = e_conf.condnames_ecaltrk_sigma;
280 
281  unsigned int encor = e_condnames_ecalonly_mean.size();
282  for (unsigned int icor=0; icor<encor; ++icor) {
283  evs.get<GBRDWrapperRcd>().get(e_condnames_ecalonly_mean[icor], forestDEH);
284  e_forestH_mean_[icor] = forestDEH.product();
285  evs.get<GBRDWrapperRcd>().get(e_condnames_ecalonly_sigma[icor], forestDEH);
286  e_forestH_sigma_[icor] = forestDEH.product();
287  }
288  for (unsigned int icor=0; icor<encor; ++icor) {
289  evs.get<GBRDWrapperRcd>().get(e_condnames_ecaltrk_mean[icor], forestDEH);
290  e_forestH_mean_[icor+encor] = forestDEH.product();
291  evs.get<GBRDWrapperRcd>().get(e_condnames_ecaltrk_sigma[icor], forestDEH);
292  e_forestH_sigma_[icor+encor] = forestDEH.product();
293  }
294 
295 }
296 
297 namespace {
298  template<typename T, typename U, typename V>
299  inline void make_consumes(T& tag,U& tok,V& sume) {
300  if(!(empty_tag == tag))
301  tok = sume.template consumes<edm::ValueMap<float> >(tag);
302  }
303 
304  template<typename T, typename U, typename V>
305  inline void make_int_consumes(T& tag,U& tok,V& sume) {
306  if(!(empty_tag == tag))
307  tok = sume.template consumes<edm::ValueMap<int> >(tag);
308  }
309 }
310 
312 
313  rhoToken_ = sumes.consumes<double>(rhoTag_);
314 
315  //setup electrons
316  if(!(empty_tag == e_conf.electron_src))
318 
319  for ( std::unordered_map<std::string, ValMapFloatTagTokenPair>::iterator imap = e_conf.tag_float_token_map.begin();
320  imap != e_conf.tag_float_token_map.end();
321  imap++) {
322  make_consumes(imap->second.first, imap->second.second, sumes);
323  }
324 
325  for ( std::unordered_map<std::string, ValMapIntTagTokenPair>::iterator imap = e_conf.tag_int_token_map.begin();
326  imap != e_conf.tag_int_token_map.end();
327  imap++) {
328  make_int_consumes(imap->second.first, imap->second.second, sumes);
329  }
330 
331  // setup photons
332  if(!(empty_tag == ph_conf.photon_src))
334 
335  for ( std::unordered_map<std::string, ValMapFloatTagTokenPair>::iterator imap = ph_conf.tag_float_token_map.begin();
336  imap != ph_conf.tag_float_token_map.end();
337  imap++) {
338  make_consumes(imap->second.first, imap->second.second, sumes);
339  }
340 
341  for ( std::unordered_map<std::string, ValMapIntTagTokenPair>::iterator imap = ph_conf.tag_int_token_map.begin();
342  imap != ph_conf.tag_int_token_map.end();
343  imap++) {
344  make_int_consumes(imap->second.first, imap->second.second, sumes);
345  }
346 }
347 
348 namespace {
349  template<typename T, typename U, typename V, typename Z>
350  inline void assignValue(const T& ptr, const U& tok, const V& map, Z& value) {
351  if( !tok.isUninitialized() ) value = map.find(tok.index())->second->get(ptr.id(),ptr.key());
352  }
353 }
354 
356 
357  // regression calculation needs no additional valuemaps
358 
359  const reco::SuperClusterRef& the_sc = ele.superCluster();
360  const edm::Ptr<reco::CaloCluster>& theseed = the_sc->seed();
361 
362  // skip HGCAL for now
363  if( theseed->seed().det() == DetId::Forward ) return;
364 
365  const int numberOfClusters = the_sc->clusters().size();
366  const bool missing_clusters = !the_sc->clusters()[numberOfClusters-1].isAvailable();
367  if( missing_clusters ) return ; // do not apply corrections in case of missing info (slimmed MiniAOD electrons)
368 
369  const bool iseb = ele.isEB();
370 
371  std::array<float, 32> eval;
372  const double raw_energy = the_sc->rawEnergy();
373  const double raw_es_energy = the_sc->preshowerEnergy();
374  const auto& full5x5_ess = ele.full5x5_showerShape();
375 
376  float e5x5Inverse = full5x5_ess.e5x5 != 0. ? vdt::fast_inv(full5x5_ess.e5x5) : 0.;
377 
378  eval[0] = raw_energy;
379  eval[1] = the_sc->etaWidth();
380  eval[2] = the_sc->phiWidth();
381  eval[3] = the_sc->seed()->energy()/raw_energy;
382  eval[4] = full5x5_ess.e5x5/raw_energy;
383  eval[5] = ele.hcalOverEcalBc();
384  eval[6] = rhoValue_;
385  eval[7] = theseed->eta() - the_sc->position().Eta();
386  eval[8] = reco::deltaPhi( theseed->phi(),the_sc->position().Phi());
387  eval[9] = full5x5_ess.r9;
388  eval[10] = full5x5_ess.sigmaIetaIeta;
389  eval[11] = full5x5_ess.sigmaIetaIphi;
390  eval[12] = full5x5_ess.sigmaIphiIphi;
391  eval[13] = full5x5_ess.eMax*e5x5Inverse;
392  eval[14] = full5x5_ess.e2nd*e5x5Inverse;
393  eval[15] = full5x5_ess.eTop*e5x5Inverse;
394  eval[16] = full5x5_ess.eBottom*e5x5Inverse;
395  eval[17] = full5x5_ess.eLeft*e5x5Inverse;
396  eval[18] = full5x5_ess.eRight*e5x5Inverse;
397  eval[19] = full5x5_ess.e2x5Max*e5x5Inverse;
398  eval[20] = full5x5_ess.e2x5Left*e5x5Inverse;
399  eval[21] = full5x5_ess.e2x5Right*e5x5Inverse;
400  eval[22] = full5x5_ess.e2x5Top*e5x5Inverse;
401  eval[23] = full5x5_ess.e2x5Bottom*e5x5Inverse;
402  eval[24] = ele.nSaturatedXtals();
403  eval[25] = std::max(0,numberOfClusters);
404 
405  // calculate coordinate variables
406  EcalClusterLocal _ecalLocal;
407  if (iseb) {
408 
409  float dummy;
410  int ieta;
411  int iphi;
412  _ecalLocal.localCoordsEB(*theseed, *iSetup_, dummy, dummy, ieta, iphi, dummy, dummy);
413  eval[26] = ieta;
414  eval[27] = iphi;
415  int signieta = ieta > 0 ? +1 : -1;
416  eval[28] = (ieta-signieta)%5;
417  eval[29] = (iphi-1)%2;
418  eval[30] = (abs(ieta)<=25)*((ieta-signieta)) + (abs(ieta)>25)*((ieta-26*signieta)%20);
419  eval[31] = (iphi-1)%20;
420 
421  } else {
422 
423  float dummy;
424  int ix;
425  int iy;
426  _ecalLocal.localCoordsEE(*theseed, *iSetup_, dummy, dummy, ix, iy, dummy, dummy);
427  eval[26] = ix;
428  eval[27] = iy;
429  eval[28] = raw_es_energy/raw_energy;
430 
431  }
432 
433  //magic numbers for MINUIT-like transformation of BDT output onto limited range
434  //(These should be stored inside the conditions object in the future as well)
435  constexpr double meanlimlow = -1.0;
436  constexpr double meanlimhigh = 3.0;
437  constexpr double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
438  constexpr double meanscale = 0.5*(meanlimhigh-meanlimlow);
439 
440  constexpr double sigmalimlow = 0.0002;
441  constexpr double sigmalimhigh = 0.5;
442  constexpr double sigmaoffset = sigmalimlow + 0.5*(sigmalimhigh-sigmalimlow);
443  constexpr double sigmascale = 0.5*(sigmalimhigh-sigmalimlow);
444 
445 
446  size_t coridx = 0;
447  float raw_pt = raw_energy*the_sc->position().rho()/the_sc->position().r();
448 
449  if (iseb && raw_pt < lowEnergy_ECALonlyThr_)
450  coridx = 0;
451  else if (iseb && raw_pt >= lowEnergy_ECALonlyThr_)
452  coridx = 1;
453  else if (!iseb && raw_pt < lowEnergy_ECALonlyThr_)
454  coridx = 2;
455  else if (!iseb && raw_pt >= lowEnergy_ECALonlyThr_)
456  coridx = 3;
457 
458  //these are the actual BDT responses
459  double rawmean = e_forestH_mean_[coridx]->GetResponse(eval.data());
460  double rawsigma = e_forestH_sigma_[coridx]->GetResponse(eval.data());
461 
462  //apply transformation to limited output range (matching the training)
463  double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
464  double sigma = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma);
465 
466  // Correct the energy. A negative energy means that the correction went
467  // outside the boundaries of the training. In this case uses raw.
468  // The resolution estimation, on the other hand should be ok.
469  if (mean < 0.) mean = 1.0;
470 
471  const double ecor = mean*(raw_energy + raw_es_energy);
472  const double sigmacor = sigma*ecor;
473 
474  ele.setCorrectedEcalEnergy(ecor);
475  ele.setCorrectedEcalEnergyError(sigmacor);
476 
477  double combinedEnergy = ecor;
478  double combinedEnergyError = sigmacor;
479 
480  auto el_track = ele.gsfTrack();
481  const float trkMomentum = el_track->pMode();
482  const float trkEta = el_track->etaMode();
483  const float trkPhi = el_track->phiMode();
484  const float trkMomentumError = std::abs(el_track->qoverpModeError())*trkMomentum*trkMomentum;
485 
486  const float eOverP = (raw_energy+raw_es_energy)*mean/trkMomentum;
487  const float fbrem = ele.fbrem();
488 
489  // E-p combination
490  if (ecor < highEnergy_ECALTRKThr_ &&
491  eOverP > eOverP_ECALTRKThr_ &&
492  std::abs(ecor - trkMomentum) < epDiffSig_ECALTRKThr_*std::sqrt(trkMomentumError*trkMomentumError+sigmacor*sigmacor) &&
493  trkMomentumError < epSig_ECALTRKThr_*trkMomentum) {
494 
495  raw_pt = ecor/cosh(trkEta);
496  if (iseb && raw_pt < lowEnergy_ECALTRKThr_)
497  coridx = 4;
498  else if (iseb && raw_pt >= lowEnergy_ECALTRKThr_)
499  coridx = 5;
500  else if (!iseb && raw_pt < lowEnergy_ECALTRKThr_)
501  coridx = 6;
502  else if (!iseb && raw_pt >= lowEnergy_ECALTRKThr_)
503  coridx = 7;
504 
505  eval[0] = ecor;
506  eval[1] = sigma/mean;
507  eval[2] = trkMomentumError/trkMomentum;
508  eval[3] = eOverP;
509  eval[4] = ele.ecalDrivenSeed();
510  eval[5] = full5x5_ess.r9;
511  eval[6] = fbrem;
512  eval[7] = trkEta;
513  eval[8] = trkPhi;
514 
515  float ecalEnergyVar = (raw_energy + raw_es_energy)*sigma;
516  float rawcombNormalization = (trkMomentumError*trkMomentumError + ecalEnergyVar*ecalEnergyVar);
517  float rawcomb = ( ecor*trkMomentumError*trkMomentumError + trkMomentum*ecalEnergyVar*ecalEnergyVar ) / rawcombNormalization;
518 
519  //these are the actual BDT responses
520  double rawmean_trk = e_forestH_mean_[coridx]->GetResponse(eval.data());
521  double rawsigma_trk = e_forestH_sigma_[coridx]->GetResponse(eval.data());
522 
523  //apply transformation to limited output range (matching the training)
524  double mean_trk = meanoffset + meanscale*vdt::fast_sin(rawmean_trk);
525  double sigma_trk = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma_trk);
526 
527  // Final correction
528  // A negative energy means that the correction went
529  // outside the boundaries of the training. In this case uses raw.
530  // The resolution estimation, on the other hand should be ok.
531  if (mean_trk < 0.) mean_trk = 1.0;
532 
533  combinedEnergy = mean_trk*rawcomb;
534  combinedEnergyError = sigma_trk*rawcomb;
535  }
536 
537  math::XYZTLorentzVector oldFourMomentum = ele.p4();
538  math::XYZTLorentzVector newFourMomentum = math::XYZTLorentzVector(oldFourMomentum.x()*combinedEnergy/oldFourMomentum.t(),
539  oldFourMomentum.y()*combinedEnergy/oldFourMomentum.t(),
540  oldFourMomentum.z()*combinedEnergy/oldFourMomentum.t(),
541  combinedEnergy);
542 
543  ele.correctMomentum(newFourMomentum, ele.trackMomentumError(), combinedEnergyError);
544 }
545 
547  modifyObject(static_cast<reco::GsfElectron&>(ele));
548 }
549 
551  // regression calculation needs no additional valuemaps
552 
553  const reco::SuperClusterRef& the_sc = pho.superCluster();
554  const edm::Ptr<reco::CaloCluster>& theseed = the_sc->seed();
555 
556  // skip HGCAL for now
557  if( theseed->seed().det() == DetId::Forward ) return;
558 
559  const int numberOfClusters = the_sc->clusters().size();
560  const bool missing_clusters = !the_sc->clusters()[numberOfClusters-1].isAvailable();
561  if( missing_clusters ) return ; // do not apply corrections in case of missing info (slimmed MiniAOD electrons)
562 
563  const bool iseb = pho.isEB();
564 
565  std::array<float, 32> eval;
566  const double raw_energy = the_sc->rawEnergy();
567  const double raw_es_energy = the_sc->preshowerEnergy();
568  const auto& full5x5_pss = pho.full5x5_showerShapeVariables();
569 
570  float e5x5Inverse = full5x5_pss.e5x5 != 0. ? vdt::fast_inv(full5x5_pss.e5x5) : 0.;
571 
572  eval[0] = raw_energy;
573  eval[1] = the_sc->etaWidth();
574  eval[2] = the_sc->phiWidth();
575  eval[3] = the_sc->seed()->energy()/raw_energy;
576  eval[4] = full5x5_pss.e5x5/raw_energy;
577  eval[5] = pho.hadronicOverEm();
578  eval[6] = rhoValue_;
579  eval[7] = theseed->eta() - the_sc->position().Eta();
580  eval[8] = reco::deltaPhi( theseed->phi(),the_sc->position().Phi());
581  eval[9] = pho.full5x5_r9();
582  eval[10] = full5x5_pss.sigmaIetaIeta;
583  eval[11] = full5x5_pss.sigmaIetaIphi;
584  eval[12] = full5x5_pss.sigmaIphiIphi;
585  eval[13] = full5x5_pss.maxEnergyXtal*e5x5Inverse;
586  eval[14] = full5x5_pss.e2nd*e5x5Inverse;
587  eval[15] = full5x5_pss.eTop*e5x5Inverse;
588  eval[16] = full5x5_pss.eBottom*e5x5Inverse;
589  eval[17] = full5x5_pss.eLeft*e5x5Inverse;
590  eval[18] = full5x5_pss.eRight*e5x5Inverse;
591  eval[19] = full5x5_pss.e2x5Max*e5x5Inverse;
592  eval[20] = full5x5_pss.e2x5Left*e5x5Inverse;
593  eval[21] = full5x5_pss.e2x5Right*e5x5Inverse;
594  eval[22] = full5x5_pss.e2x5Top*e5x5Inverse;
595  eval[23] = full5x5_pss.e2x5Bottom*e5x5Inverse;
596  eval[24] = pho.nSaturatedXtals();
597  eval[25] = std::max(0,numberOfClusters);
598 
599  // calculate coordinate variables
600  EcalClusterLocal _ecalLocal;
601 
602  if (iseb) {
603 
604  float dummy;
605  int ieta;
606  int iphi;
607  _ecalLocal.localCoordsEB(*theseed, *iSetup_, dummy, dummy, ieta, iphi, dummy, dummy);
608  eval[26] = ieta;
609  eval[27] = iphi;
610  int signieta = ieta > 0 ? +1 : -1;
611  eval[28] = (ieta-signieta)%5;
612  eval[29] = (iphi-1)%2;
613  eval[30] = (abs(ieta)<=25)*((ieta-signieta)) + (abs(ieta)>25)*((ieta-26*signieta)%20);
614  eval[31] = (iphi-1)%20;
615 
616  } else {
617 
618  float dummy;
619  int ix;
620  int iy;
621  _ecalLocal.localCoordsEE(*theseed, *iSetup_, dummy, dummy, ix, iy, dummy, dummy);
622  eval[26] = ix;
623  eval[27] = iy;
624  eval[28] = raw_es_energy/raw_energy;
625 
626  }
627 
628  //magic numbers for MINUIT-like transformation of BDT output onto limited range
629  //(These should be stored inside the conditions object in the future as well)
630  constexpr double meanlimlow = -1.0;
631  constexpr double meanlimhigh = 3.0;
632  constexpr double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
633  constexpr double meanscale = 0.5*(meanlimhigh-meanlimlow);
634 
635  constexpr double sigmalimlow = 0.0002;
636  constexpr double sigmalimhigh = 0.5;
637  constexpr double sigmaoffset = sigmalimlow + 0.5*(sigmalimhigh-sigmalimlow);
638  constexpr double sigmascale = 0.5*(sigmalimhigh-sigmalimlow);
639 
640  size_t coridx = 0;
641  float raw_pt = raw_energy*the_sc->position().rho()/the_sc->position().r();
642 
643  if (iseb && raw_pt < lowEnergy_ECALonlyThr_)
644  coridx = 0;
645  else if (iseb && raw_pt >= lowEnergy_ECALonlyThr_)
646  coridx = 1;
647  else if (!iseb && raw_pt < lowEnergy_ECALonlyThr_)
648  coridx = 2;
649  else if (!iseb && raw_pt >= lowEnergy_ECALonlyThr_)
650  coridx = 3;
651 
652  //these are the actual BDT responses
653  double rawmean = ph_forestH_mean_[coridx]->GetResponse(eval.data());
654  double rawsigma = ph_forestH_sigma_[coridx]->GetResponse(eval.data());
655 
656  //apply transformation to limited output range (matching the training)
657  double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
658  double sigma = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma);
659 
660  // Correct the energy. A negative energy means that the correction went
661  // outside the boundaries of the training. In this case uses raw.
662  // The resolution estimation, on the other hand should be ok.
663  if (mean < 0.) mean = 1.0;
664 
665  const double ecor = mean*(raw_energy + raw_es_energy);
666  const double sigmacor = sigma*ecor;
667 
668  pho.setCorrectedEnergy(reco::Photon::P4type::regression2, ecor, sigmacor, true);
669 }
670 
672  modifyObject(static_cast<reco::Photon&>(pho));
673 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool isAvailable() const
Definition: Ref.h:577
T getParameter(std::string const &) const
Analysis-level Photon class.
Definition: Photon.h:47
float nSaturatedXtals() const
Definition: GsfElectron.h:512
float trackMomentumError() const
Definition: GsfElectron.h:825
Definition: Photon.py:1
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
std::unordered_map< std::string, ValMapIntTagTokenPair > tag_int_token_map
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:225
reco::SuperClusterRef superCluster() const
Ref to SuperCluster.
key_type key() const
Definition: Ptr.h:185
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
std::unordered_map< std::string, ValMapIntTagTokenPair > tag_int_token_map
std::unordered_map< unsigned, edm::Ptr< reco::Photon > > phos_by_oop
void correctMomentum(const LorentzVector &p4, float trackMomentumError, float p4Error)
Definition: GsfElectron.h:847
void modifyObject(reco::GsfElectron &) const override final
std::unordered_map< std::string, ValMapFloatTagTokenPair > tag_float_token_map
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< const GBRForestD * > ph_forestH_sigma_
float fbrem() const
Definition: GsfElectron.h:750
std::unordered_map< std::string, ValMapFloatTagTokenPair > tag_float_token_map
std::unordered_map< unsigned, edm::Handle< edm::ValueMap< int > > > pho_int_vmaps
#define constexpr
bool isEB() const
Definition: GsfElectron.h:352
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:166
Definition: HeavyIon.h:7
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
float full5x5_r9() const
Definition: Photon.h:240
std::unordered_map< unsigned, edm::Handle< edm::ValueMap< float > > > pho_vmaps
void setCorrectedEcalEnergyError(float newEnergyError)
Definition: GsfElectron.cc:179
const std::string & name() const
std::pair< edm::InputTag, ValMapFloatToken > ValMapFloatTagTokenPair
void setEventContent(const edm::EventSetup &) override final
T sqrt(T t)
Definition: SSEVec.h:18
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:184
std::unordered_map< unsigned, edm::Ptr< reco::GsfElectron > > eles_by_oop
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: value.py:1
float hadronicOverEm() const
the total hadronic over electromagnetic fraction
Definition: Photon.h:204
std::vector< std::string > getParameterNames() const
std::unordered_map< unsigned, edm::Handle< edm::ValueMap< float > > > ele_vmaps
edm::EDGetTokenT< edm::ValueMap< int > > ValMapIntToken
const edm::Ptr< reco::Candidate > & originalObjectRef() const
reference to original object. Returns a null reference if not available
Definition: PATObject.h:500
float hcalOverEcalBc() const
Definition: GsfElectron.h:446
std::vector< const GBRForestD * > ph_forestH_mean_
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
edm::EDGetTokenT< edm::View< pat::Photon > > tok_photon_src
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:205
std::pair< edm::InputTag, ValMapIntToken > ValMapIntTagTokenPair
const T & get() const
Definition: EventSetup.h:55
edm::EDGetTokenT< edm::ValueMap< float > > ValMapFloatToken
Analysis-level electron class.
Definition: Electron.h:52
std::vector< const GBRForestD * > e_forestH_mean_
void localCoordsEB(const reco::CaloCluster &bclus, const edm::EventSetup &es, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt) const
return(e1-e2)*(e1-e2)+dp *dp
bool isEB() const
Definition: Photon.h:121
void setCorrectedEcalEnergy(float newEnergy)
Definition: GsfElectron.cc:182
const ShowerShape & full5x5_showerShape() const
Definition: GsfElectron.h:475
fixed size matrix
HLT enums.
const ShowerShape & full5x5_showerShapeVariables() const
Definition: Photon.h:198
void setEvent(const edm::Event &) override final
bool isUninitialized() const
Definition: EDGetToken.h:73
std::unordered_map< unsigned, edm::Handle< edm::ValueMap< int > > > ele_int_vmaps
#define DEFINE_EDM_PLUGIN(factory, type, name)
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:169
void localCoordsEE(const reco::CaloCluster &bclus, const edm::EventSetup &es, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
void setConsumes(edm::ConsumesCollector &) override final
long double T
T const * product() const
Definition: ESHandle.h:86
std::vector< const GBRForestD * > e_forestH_sigma_
edm::EDGetTokenT< double > rhoToken_
EGExtraInfoModifierFromDB(const edm::ParameterSet &conf)
edm::EDGetTokenT< edm::View< pat::Electron > > tok_electron_src
float nSaturatedXtals() const
Definition: Photon.h:254
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:185
bool ecalDrivenSeed() const
Definition: GsfElectron.h:188