CMS 3D CMS Logo

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