CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Attributes
EGRegressionModifierV1 Class Reference
Inheritance diagram for EGRegressionModifierV1:
ModifyObjectValueBase

Classes

struct  CondNames
 

Public Member Functions

 EGRegressionModifierV1 (const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
 
void modifyObject (reco::GsfElectron &) const final
 
void modifyObject (reco::Photon &) const final
 
void modifyObject (pat::Electron &ele) const final
 
void modifyObject (pat::Photon &pho) const final
 
void setEvent (const edm::Event &) final
 
void setEventContent (const edm::EventSetup &) final
 
- Public Member Functions inherited from ModifyObjectValueBase
virtual void modifyObject (reco::Muon &) const
 
virtual void modifyObject (reco::BaseTau &) const
 
virtual void modifyObject (reco::Jet &) const
 
virtual void modifyObject (pat::Muon &) const
 
virtual void modifyObject (pat::Tau &) const
 
virtual void modifyObject (pat::Jet &) const
 
 ModifyObjectValueBase (const edm::ParameterSet &conf)
 
const std::string & name () const
 
virtual ~ModifyObjectValueBase ()
 

Private Attributes

const bool applyExtraHighEnergyProtection_
 
const bool autoDetectBunchSpacing_
 
int bunchspacing_
 
edm::EDGetTokenT< unsigned int > bunchSpacingToken_
 
edm::ESHandle< CaloGeometrycaloGeomH_
 
std::string condNamesWeight25ns_
 
std::string condNamesWeight50ns_
 
CondNames eleCondNames_
 
std::vector< const GBRForestD * > eleForestsMean_
 
std::vector< const GBRForestD * > eleForestsSigma_
 
const GBRForestepForest_
 
int nVtx_
 
CondNames phoCondNames_
 
std::vector< const GBRForestD * > phoForestsMean_
 
std::vector< const GBRForestD * > phoForestsSigma_
 
edm::EDGetTokenT< double > rhoToken_
 
float rhoValue_
 
edm::Handle< reco::VertexCollectionvtxH_
 
edm::EDGetTokenT< reco::VertexCollectionvtxToken_
 

Detailed Description

Definition at line 18 of file EGRegressionModifierV1.cc.

Constructor & Destructor Documentation

EGRegressionModifierV1::EGRegressionModifierV1 ( const edm::ParameterSet conf,
edm::ConsumesCollector cc 
)

Definition at line 68 of file EGRegressionModifierV1.cc.

References autoDetectBunchSpacing_, bunchSpacingToken_, condNamesWeight25ns_, condNamesWeight50ns_, edm::ConsumesCollector::consumes(), eleCondNames_, lowPtElectrons_cff::electrons, edm::ParameterSet::getParameter(), EGRegressionModifierV1::CondNames::mean50ns, phoCondNames_, muons_cff::photons, and AlCaHLTBitMon_QueryRunRegistry::string.

69  : ModifyObjectValueBase(conf)
70  , autoDetectBunchSpacing_(conf.getParameter<bool>("autoDetectBunchSpacing"))
71  , bunchspacing_(autoDetectBunchSpacing_ ? 450 : conf.getParameter<int>("manualBunchSpacing"))
72  , rhoToken_(cc.consumes<double>(conf.getParameter<edm::InputTag>("rhoCollection")))
73  , vtxToken_(cc.consumes<reco::VertexCollection>(conf.getParameter<edm::InputTag>("vertexCollection")))
74  , applyExtraHighEnergyProtection_(conf.getParameter<bool>("applyExtraHighEnergyProtection"))
75 {
77  bunchSpacingToken_ = cc.consumes<unsigned int>(conf.getParameter<edm::InputTag>("bunchSpacingTag"));
78 
79  const edm::ParameterSet& electrons = conf.getParameter<edm::ParameterSet>("electron_config");
80  eleCondNames_ = CondNames {
81  .mean50ns = electrons.getParameter<std::vector<std::string> >("regressionKey_50ns"),
82  .sigma50ns = electrons.getParameter<std::vector<std::string> >("uncertaintyKey_50ns"),
83  .mean25ns = electrons.getParameter<std::vector<std::string> >("regressionKey_25ns"),
84  .sigma25ns = electrons.getParameter<std::vector<std::string> >("uncertaintyKey_25ns")
85  };
86  condNamesWeight50ns_ = electrons.getParameter<std::string>("combinationKey_50ns");
87  condNamesWeight25ns_ = electrons.getParameter<std::string>("combinationKey_25ns");
88 
89  const edm::ParameterSet& photons = conf.getParameter<edm::ParameterSet>("photon_config");
90  phoCondNames_ = CondNames {
91  .mean50ns = photons.getParameter<std::vector<std::string>>("regressionKey_50ns"),
92  .sigma50ns = photons.getParameter<std::vector<std::string>>("uncertaintyKey_50ns"),
93  .mean25ns = photons.getParameter<std::vector<std::string>>("regressionKey_25ns"),
94  .sigma25ns = photons.getParameter<std::vector<std::string>>("uncertaintyKey_25ns")
95  };
96 }
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< double > rhoToken_
ModifyObjectValueBase(const edm::ParameterSet &conf)
edm::EDGetTokenT< unsigned int > bunchSpacingToken_

Member Function Documentation

void EGRegressionModifierV1::modifyObject ( reco::GsfElectron ele) const
finalvirtual

Reimplemented from ModifyObjectValueBase.

Definition at line 130 of file EGRegressionModifierV1.cc.

References funct::abs(), applyExtraHighEnergyProtection_, caloGeomH_, reco::GsfElectron::classification(), constexpr, reco::GsfElectron::correctMomentum(), reco::deltaPhi(), reco::deltaR(), reco::GsfElectron::ecalDriven(), eleForestsMean_, eleForestsSigma_, reco::CaloCluster::energy(), epForest_, reco::CaloCluster::eta(), GBRForest::GetResponse(), createfilelist::int, edm::Ref< C, T, F >::isAvailable(), reco::GsfElectron::isEB(), egammaTools::localEcalClusterCoordsEB(), egammaTools::localEcalClusterCoordsEE(), SiStripPI::max, jets_cff::maxDR, SiStripPI::mean, nVtx_, reco::GsfElectron::p4(), reco::CaloCluster::phi(), mathSSE::return(), reco::GsfElectron::setCorrectedEcalEnergy(), reco::GsfElectron::setCorrectedEcalEnergyError(), reco::GsfElectron::showerShape(), mathSSE::sqrt(), reco::GsfElectron::superCluster(), reco::GsfElectron::trackerDrivenSeed(), reco::GsfElectron::trackMomentumAtVtx(), reco::GsfElectron::trackMomentumError(), and mps_merge::weight.

130  {
131  // regression calculation needs no additional valuemaps
132 
133  const reco::SuperClusterRef& superClus = ele.superCluster();
134  const edm::Ptr<reco::CaloCluster>& theseed = superClus->seed();
135  const int numberOfClusters = superClus->clusters().size();
136  const bool missing_clusters = !superClus->clusters()[numberOfClusters-1].isAvailable();
137 
138  if( missing_clusters ) return ; // do not apply corrections in case of missing info (slimmed MiniAOD electrons)
139 
140  std::array<float, 33> eval;
141  const double rawEnergy = superClus->rawEnergy();
142  const auto& ess = ele.showerShape();
143 
144  // SET INPUTS
145  eval[0] = nVtx_;
146  eval[1] = rawEnergy;
147  eval[2] = superClus->eta();
148  eval[3] = superClus->phi();
149  eval[4] = superClus->etaWidth();
150  eval[5] = superClus->phiWidth();
151  eval[6] = ess.r9;
152  eval[7] = theseed->energy()/rawEnergy;
153  eval[8] = ess.eMax/rawEnergy;
154  eval[9] = ess.e2nd/rawEnergy;
155  eval[10] = (ess.eLeft + ess.eRight != 0.f ? (ess.eLeft-ess.eRight)/(ess.eLeft+ess.eRight) : 0.f);
156  eval[11] = (ess.eTop + ess.eBottom != 0.f ? (ess.eTop-ess.eBottom)/(ess.eTop+ess.eBottom) : 0.f);
157  eval[12] = ess.sigmaIetaIeta;
158  eval[13] = ess.sigmaIetaIphi;
159  eval[14] = ess.sigmaIphiIphi;
160  eval[15] = std::max(0,numberOfClusters-1);
161 
162  // calculate sub-cluster variables
163  std::vector<float> clusterRawEnergy;
164  clusterRawEnergy.resize(std::max(3, numberOfClusters), 0);
165  std::vector<float> clusterDEtaToSeed;
166  clusterDEtaToSeed.resize(std::max(3, numberOfClusters), 0);
167  std::vector<float> clusterDPhiToSeed;
168  clusterDPhiToSeed.resize(std::max(3, numberOfClusters), 0);
169  float clusterMaxDR = 999.;
170  float clusterMaxDRDPhi = 999.;
171  float clusterMaxDRDEta = 999.;
172  float clusterMaxDRRawEnergy = 0.;
173 
174  size_t iclus = 0;
175  float maxDR = 0;
176  // loop over all clusters that aren't the seed
177  for (auto const& pclus : superClus->clusters())
178  {
179  if(theseed == pclus )
180  continue;
181  clusterRawEnergy[iclus] = pclus->energy();
182  clusterDPhiToSeed[iclus] = reco::deltaPhi(pclus->phi(),theseed->phi());
183  clusterDEtaToSeed[iclus] = pclus->eta() - theseed->eta();
184 
185  // find cluster with max dR
186  const auto the_dr = reco::deltaR(*pclus, *theseed);
187  if(the_dr > maxDR) {
188  maxDR = the_dr;
189  clusterMaxDR = maxDR;
190  clusterMaxDRDPhi = clusterDPhiToSeed[iclus];
191  clusterMaxDRDEta = clusterDEtaToSeed[iclus];
192  clusterMaxDRRawEnergy = clusterRawEnergy[iclus];
193  }
194  ++iclus;
195  }
196 
197  eval[16] = clusterMaxDR;
198  eval[17] = clusterMaxDRDPhi;
199  eval[18] = clusterMaxDRDEta;
200  eval[19] = clusterMaxDRRawEnergy/rawEnergy;
201  eval[20] = clusterRawEnergy[0]/rawEnergy;
202  eval[21] = clusterRawEnergy[1]/rawEnergy;
203  eval[22] = clusterRawEnergy[2]/rawEnergy;
204  eval[23] = clusterDPhiToSeed[0];
205  eval[24] = clusterDPhiToSeed[1];
206  eval[25] = clusterDPhiToSeed[2];
207  eval[26] = clusterDEtaToSeed[0];
208  eval[27] = clusterDEtaToSeed[1];
209  eval[28] = clusterDEtaToSeed[2];
210 
211  // calculate coordinate variables
212  const bool isEB = ele.isEB();
213  float dummy;
214  int iPhi;
215  int iEta;
216  float cryPhi;
217  float cryEta;
218  if (ele.isEB())
219  egammaTools::localEcalClusterCoordsEB(*theseed, *caloGeomH_, cryEta, cryPhi, iEta, iPhi, dummy, dummy);
220  else
221  egammaTools::localEcalClusterCoordsEE(*theseed, *caloGeomH_, cryEta, cryPhi, iEta, iPhi, dummy, dummy);
222 
223  if (isEB) {
224  eval[29] = cryEta;
225  eval[30] = cryPhi;
226  eval[31] = iEta;
227  eval[32] = iPhi;
228  } else {
229  eval[29] = superClus->preshowerEnergy()/superClus->rawEnergy();
230  }
231 
232  //magic numbers for MINUIT-like transformation of BDT output onto limited range
233  //(These should be stored inside the conditions object in the future as well)
234  constexpr double meanlimlow = 0.2;
235  constexpr double meanlimhigh = 2.0;
236  constexpr double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
237  constexpr double meanscale = 0.5*(meanlimhigh-meanlimlow);
238 
239  constexpr double sigmalimlow = 0.0002;
240  constexpr double sigmalimhigh = 0.5;
241  constexpr double sigmaoffset = sigmalimlow + 0.5*(sigmalimhigh-sigmalimlow);
242  constexpr double sigmascale = 0.5*(sigmalimhigh-sigmalimlow);
243 
244  const int coridx = isEB ? 0 : 1;
245 
246  //these are the actual BDT responses
247  double rawmean = eleForestsMean_[coridx]->GetResponse(eval.data());
248  double rawsigma = eleForestsSigma_[coridx]->GetResponse(eval.data());
249 
250  //apply transformation to limited output range (matching the training)
251  double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
252  double sigma = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma);
253 
254  //regression target is ln(Etrue/Eraw)
255  //so corrected energy is ecor=exp(mean)*e, uncertainty is exp(mean)*eraw*sigma=ecor*sigma
256  double ecor = mean*(eval[1]);
257  if (!isEB)
258  ecor = mean*(eval[1]+superClus->preshowerEnergy());
259  const double sigmacor = sigma*ecor;
260 
261  ele.setCorrectedEcalEnergy(ecor);
262  ele.setCorrectedEcalEnergyError(sigmacor);
263 
264  // E-p combination
265  std::array<float, 11> eval_ep;
266 
267  const float ep = ele.trackMomentumAtVtx().R();
268  const float tot_energy = superClus->rawEnergy()+superClus->preshowerEnergy();
269  const float momentumError = ele.trackMomentumError();
270  const float trkMomentumRelError = ele.trackMomentumError()/ep;
271  const float eOverP = tot_energy*mean/ep;
272  eval_ep[0] = tot_energy*mean;
273  eval_ep[1] = sigma/mean;
274  eval_ep[2] = ep;
275  eval_ep[3] = trkMomentumRelError;
276  eval_ep[4] = sigma/mean/trkMomentumRelError;
277  eval_ep[5] = tot_energy*mean/ep;
278  eval_ep[6] = tot_energy*mean/ep*sqrt(sigma/mean*sigma/mean+trkMomentumRelError*trkMomentumRelError);
279  eval_ep[7] = ele.ecalDriven();
280  eval_ep[8] = ele.trackerDrivenSeed();
281  eval_ep[9] = int(ele.classification());//eleClass;
282  eval_ep[10] = isEB;
283 
284  // CODE FOR FUTURE SEMI_PARAMETRIC
285  //double rawweight = ep_forestsMean_[coridx]->GetResponse(eval_ep.data());
287  //double weight = meanoffset + meanscale*vdt::fast_sin(rawweight);
289 
290  // CODE FOR STANDARD BDT
291  double weight = 0.;
292  if ( eOverP > 0.025 &&
293  std::abs(ep-ecor) < 15.*std::sqrt( momentumError*momentumError + sigmacor*sigmacor ) &&
294  (!applyExtraHighEnergyProtection_ || ((momentumError < 10.*ep) || (ecor < 200.)))
295  ) {
296  // protection against crazy track measurement
297  weight = std::clamp(epForest_->GetResponse(eval_ep.data()), 0., 1.);
298  }
299 
300  double combinedMomentum = weight*ele.trackMomentumAtVtx().R() + (1.-weight)*ecor;
301  double combinedMomentumError = sqrt(weight*weight*ele.trackMomentumError()*ele.trackMomentumError() + (1.-weight)*(1.-weight)*sigmacor*sigmacor);
302 
303  math::XYZTLorentzVector oldMomentum = ele.p4();
304  math::XYZTLorentzVector newMomentum( oldMomentum.x()*combinedMomentum/oldMomentum.t(),
305  oldMomentum.y()*combinedMomentum/oldMomentum.t(),
306  oldMomentum.z()*combinedMomentum/oldMomentum.t(),
307  combinedMomentum );
308 
309  ele.correctMomentum(newMomentum, ele.trackMomentumError(), combinedMomentumError);
310 }
double GetResponse(const float *vector) const
Definition: GBRForest.h:53
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
const ShowerShape & showerShape() const
Definition: GsfElectron.h:458
bool isAvailable() const
Definition: Ref.h:575
void localEcalClusterCoordsEB(const reco::CaloCluster &bclus, const CaloGeometry &geom, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt)
float trackMomentumError() const
Definition: GsfElectron.h:848
const LorentzVector & p4(P4Kind kind) const
Definition: GsfElectron.cc:228
edm::ESHandle< CaloGeometry > caloGeomH_
void correctMomentum(const LorentzVector &p4, float trackMomentumError, float p4Error)
Definition: GsfElectron.h:870
math::XYZVectorF trackMomentumAtVtx() const
Definition: GsfElectron.h:295
Definition: weight.py:1
bool isEB() const
Definition: GsfElectron.h:356
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:168
return((rh^lh)&mask)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
void setCorrectedEcalEnergyError(float newEnergyError)
Definition: GsfElectron.cc:179
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< const GBRForestD * > eleForestsSigma_
std::vector< const GBRForestD * > eleForestsMean_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:126
void localEcalClusterCoordsEE(const reco::CaloCluster &bclus, const CaloGeometry &geom, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt)
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
Classification classification() const
Definition: GsfElectron.h:768
void setCorrectedEcalEnergy(float newEnergy)
Definition: GsfElectron.cc:182
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:185
bool trackerDrivenSeed() const
Definition: GsfElectron.h:189
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:171
bool ecalDriven() const
Definition: GsfElectron.cc:174
#define constexpr
void EGRegressionModifierV1::modifyObject ( reco::Photon pho) const
finalvirtual

this is 1*abs(ieta)/ieta in original training

duplicated variables but this was trained like that

duplicated variables but this was trained like that

Reimplemented from ModifyObjectValueBase.

Definition at line 312 of file EGRegressionModifierV1.cc.

References funct::abs(), reco::deltaPhi(), reco::Photon::e5x5(), reco::CaloCluster::energy(), reco::CaloCluster::eta(), reco::Photon::hadronicOverEm(), edm::Ref< C, T, F >::isAvailable(), reco::Photon::isEB(), SiStripPI::max, SiStripPI::mean, nVtx_, reco::CaloCluster::phi(), phoForestsMean_, phoForestsSigma_, reco::Photon::r9(), mathSSE::return(), rhoValue_, reco::CaloCluster::seed(), reco::Photon::setCorrectedEnergy(), reco::Photon::showerShapeVariables(), and reco::Photon::superCluster().

312  {
313  // regression calculation needs no additional valuemaps
314 
315  std::array<float, 35> eval;
316  const reco::SuperClusterRef& superClus = pho.superCluster();
317  const edm::Ptr<reco::CaloCluster>& theseed = superClus->seed();
318 
319  const int numberOfClusters = superClus->clusters().size();
320  const bool missing_clusters = !superClus->clusters()[numberOfClusters-1].isAvailable();
321 
322  if( missing_clusters ) return ; // do not apply corrections in case of missing info (slimmed MiniAOD electrons)
323 
324  const double rawEnergy = superClus->rawEnergy();
325  const auto& ess = pho.showerShapeVariables();
326 
327  // SET INPUTS
328  eval[0] = rawEnergy;
329  eval[1] = pho.r9();
330  eval[2] = superClus->etaWidth();
331  eval[3] = superClus->phiWidth();
332  eval[4] = std::max(0,numberOfClusters - 1);
333  eval[5] = pho.hadronicOverEm();
334  eval[6] = rhoValue_;
335  eval[7] = nVtx_;
336  eval[8] = theseed->eta()-superClus->position().Eta();
337  eval[9] = reco::deltaPhi(theseed->phi(),superClus->position().Phi());
338  eval[10] = theseed->energy()/rawEnergy;
339  eval[11] = ess.e3x3/ess.e5x5;
340  eval[12] = ess.sigmaIetaIeta;
341  eval[13] = ess.sigmaIphiIphi;
342  eval[14] = ess.sigmaIetaIphi/(ess.sigmaIphiIphi*ess.sigmaIetaIeta);
343  eval[15] = ess.maxEnergyXtal/ess.e5x5;
344  eval[16] = ess.e2nd/ess.e5x5;
345  eval[17] = ess.eTop/ess.e5x5;
346  eval[18] = ess.eBottom/ess.e5x5;
347  eval[19] = ess.eLeft/ess.e5x5;
348  eval[20] = ess.eRight/ess.e5x5;
349  eval[21] = ess.e2x5Max/ess.e5x5;
350  eval[22] = ess.e2x5Left/ess.e5x5;
351  eval[23] = ess.e2x5Right/ess.e5x5;
352  eval[24] = ess.e2x5Top/ess.e5x5;
353  eval[25] = ess.e2x5Bottom/ess.e5x5;
354 
355  const bool isEB = pho.isEB();
356  if (isEB) {
357  EBDetId ebseedid(theseed->seed());
358  eval[26] = pho.e5x5()/theseed->energy();
359  int ieta = ebseedid.ieta();
360  int iphi = ebseedid.iphi();
361  eval[27] = ieta;
362  eval[28] = iphi;
363  int signieta = ieta > 0 ? +1 : -1;
364  eval[29] = (ieta-signieta)%5;
365  eval[30] = (iphi-1)%2;
366  eval[31] = (abs(ieta)<=25)*((ieta-signieta)) + (abs(ieta)>25)*((ieta-26*signieta)%20);
367  eval[32] = (iphi-1)%20;
368  eval[33] = ieta;
369  eval[34] = iphi;
370  } else {
371  EEDetId eeseedid(theseed->seed());
372  eval[26] = superClus->preshowerEnergy()/rawEnergy;
373  eval[27] = superClus->preshowerEnergyPlane1()/rawEnergy;
374  eval[28] = superClus->preshowerEnergyPlane2()/rawEnergy;
375  eval[29] = eeseedid.ix();
376  eval[30] = eeseedid.iy();
377  }
378 
379  //magic numbers for MINUIT-like transformation of BDT output onto limited range
380  //(These should be stored inside the conditions object in the future as well)
381  const double meanlimlow = 0.2;
382  const double meanlimhigh = 2.0;
383  const double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
384  const double meanscale = 0.5*(meanlimhigh-meanlimlow);
385 
386  const double sigmalimlow = 0.0002;
387  const double sigmalimhigh = 0.5;
388  const double sigmaoffset = sigmalimlow + 0.5*(sigmalimhigh-sigmalimlow);
389  const double sigmascale = 0.5*(sigmalimhigh-sigmalimlow);
390 
391  const int coridx = isEB ? 0 : 1;
392 
393  //these are the actual BDT responses
394  const double rawmean = phoForestsMean_[coridx]->GetResponse(eval.data());
395  const double rawsigma = phoForestsSigma_[coridx]->GetResponse(eval.data());
396  //apply transformation to limited output range (matching the training)
397  const double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
398  const double sigma = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma);
399 
400  //regression target is ln(Etrue/Eraw)
401  //so corrected energy is ecor=exp(mean)*e, uncertainty is exp(mean)*eraw*sigma=ecor*sigma
402  const double ecor = isEB ? mean*eval[0] : mean*(eval[0]+superClus->preshowerEnergy());
403 
404  const double sigmacor = sigma*ecor;
405  pho.setCorrectedEnergy(reco::Photon::P4type::regression2, ecor, sigmacor, true);
406 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
bool isAvailable() const
Definition: Ref.h:575
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
float e5x5() const
Definition: Photon.h:234
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:168
return((rh^lh)&mask)
std::vector< const GBRForestD * > phoForestsMean_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:126
float hadronicOverEm() const
the total hadronic over electromagnetic fraction
Definition: Photon.h:212
std::vector< const GBRForestD * > phoForestsSigma_
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:207
bool isEB() const
Definition: Photon.h:121
float r9() const
Definition: Photon.h:240
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:171
const ShowerShape & showerShapeVariables() const
Definition: Photon.h:205
void EGRegressionModifierV1::modifyObject ( pat::Electron ele) const
inlinefinalvirtual

Reimplemented from ModifyObjectValueBase.

Definition at line 30 of file EGRegressionModifierV1.cc.

References modifyObject().

Referenced by modifyObject().

30 { modifyObject(static_cast<reco::GsfElectron&>(ele)); }
void modifyObject(reco::GsfElectron &) const final
void EGRegressionModifierV1::modifyObject ( pat::Photon pho) const
inlinefinalvirtual

Reimplemented from ModifyObjectValueBase.

Definition at line 31 of file EGRegressionModifierV1.cc.

References modifyObject().

Referenced by modifyObject().

31 { modifyObject(static_cast<reco::Photon&>(pho)); }
void modifyObject(reco::GsfElectron &) const final
void EGRegressionModifierV1::setEvent ( const edm::Event evt)
finalvirtual

Reimplemented from ModifyObjectValueBase.

Definition at line 98 of file EGRegressionModifierV1.cc.

References autoDetectBunchSpacing_, bunchspacing_, bunchSpacingToken_, edm::Event::getByToken(), nVtx_, rhoToken_, rhoValue_, vtxH_, and vtxToken_.

99 {
101  edm::Handle<unsigned int> bunchSpacingH;
102  evt.getByToken(bunchSpacingToken_,bunchSpacingH);
103  bunchspacing_ = *bunchSpacingH;
104  }
105 
106  edm::Handle<double> rhoH;
107  evt.getByToken(rhoToken_, rhoH);
108  rhoValue_ = *rhoH;
109 
110  evt.getByToken(vtxToken_, vtxH_);
111  nVtx_ = vtxH_->size();
112 }
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
edm::Handle< reco::VertexCollection > vtxH_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< double > rhoToken_
edm::EDGetTokenT< unsigned int > bunchSpacingToken_
void EGRegressionModifierV1::setEventContent ( const edm::EventSetup evs)
finalvirtual

Reimplemented from ModifyObjectValueBase.

Definition at line 114 of file EGRegressionModifierV1.cc.

References bunchspacing_, caloGeomH_, condNamesWeight25ns_, condNamesWeight50ns_, eleCondNames_, eleForestsMean_, eleForestsSigma_, epForest_, edm::EventSetup::get(), EGRegressionModifierV1::CondNames::mean25ns, EGRegressionModifierV1::CondNames::mean50ns, phoCondNames_, phoForestsMean_, phoForestsSigma_, edm::ESHandle< T >::product(), retrieveGBRForests(), EGRegressionModifierV1::CondNames::sigma25ns, EGRegressionModifierV1::CondNames::sigma50ns, and AlCaHLTBitMon_QueryRunRegistry::string.

115 {
116  evs.get<CaloGeometryRecord>().get(caloGeomH_);
117 
120 
123 
124  edm::ESHandle<GBRForest> forestEH;
125  const std::string ep_condnames_weight = (bunchspacing_ == 25) ? condNamesWeight25ns_ : condNamesWeight50ns_;
126  evs.get<GBRWrapperRcd>().get(ep_condnames_weight, forestEH);
127  epForest_ = forestEH.product();
128 }
std::vector< std::string > sigma25ns
edm::ESHandle< CaloGeometry > caloGeomH_
std::vector< const GBRForestD * > phoForestsMean_
std::vector< std::string > sigma50ns
std::vector< const GBRForestD * > retrieveGBRForests(edm::EventSetup const &evs, std::vector< std::string > const &names)
std::vector< const GBRForestD * > eleForestsSigma_
std::vector< const GBRForestD * > eleForestsMean_
std::vector< const GBRForestD * > phoForestsSigma_
T get() const
Definition: EventSetup.h:71
T const * product() const
Definition: ESHandle.h:86

Member Data Documentation

const bool EGRegressionModifierV1::applyExtraHighEnergyProtection_
private

Definition at line 57 of file EGRegressionModifierV1.cc.

Referenced by modifyObject().

const bool EGRegressionModifierV1::autoDetectBunchSpacing_
private

Definition at line 48 of file EGRegressionModifierV1.cc.

Referenced by EGRegressionModifierV1(), and setEvent().

int EGRegressionModifierV1::bunchspacing_
private

Definition at line 49 of file EGRegressionModifierV1.cc.

Referenced by setEvent(), and setEventContent().

edm::EDGetTokenT<unsigned int> EGRegressionModifierV1::bunchSpacingToken_
private

Definition at line 50 of file EGRegressionModifierV1.cc.

Referenced by EGRegressionModifierV1(), and setEvent().

edm::ESHandle<CaloGeometry> EGRegressionModifierV1::caloGeomH_
private

Definition at line 56 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

std::string EGRegressionModifierV1::condNamesWeight25ns_
private

Definition at line 46 of file EGRegressionModifierV1.cc.

Referenced by EGRegressionModifierV1(), and setEventContent().

std::string EGRegressionModifierV1::condNamesWeight50ns_
private

Definition at line 45 of file EGRegressionModifierV1.cc.

Referenced by EGRegressionModifierV1(), and setEventContent().

CondNames EGRegressionModifierV1::eleCondNames_
private

Definition at line 42 of file EGRegressionModifierV1.cc.

Referenced by EGRegressionModifierV1(), and setEventContent().

std::vector<const GBRForestD*> EGRegressionModifierV1::eleForestsMean_
private

Definition at line 61 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

std::vector<const GBRForestD*> EGRegressionModifierV1::eleForestsSigma_
private

Definition at line 62 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

const GBRForest* EGRegressionModifierV1::epForest_
private

Definition at line 63 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

int EGRegressionModifierV1::nVtx_
private

Definition at line 53 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEvent().

CondNames EGRegressionModifierV1::phoCondNames_
private

Definition at line 43 of file EGRegressionModifierV1.cc.

Referenced by EGRegressionModifierV1(), and setEventContent().

std::vector<const GBRForestD*> EGRegressionModifierV1::phoForestsMean_
private

Definition at line 59 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

std::vector<const GBRForestD*> EGRegressionModifierV1::phoForestsSigma_
private

Definition at line 60 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

edm::EDGetTokenT<double> EGRegressionModifierV1::rhoToken_
private

Definition at line 52 of file EGRegressionModifierV1.cc.

Referenced by setEvent().

float EGRegressionModifierV1::rhoValue_
private

Definition at line 51 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEvent().

edm::Handle<reco::VertexCollection> EGRegressionModifierV1::vtxH_
private

Definition at line 55 of file EGRegressionModifierV1.cc.

Referenced by setEvent().

edm::EDGetTokenT<reco::VertexCollection> EGRegressionModifierV1::vtxToken_
private

Definition at line 54 of file EGRegressionModifierV1.cc.

Referenced by setEvent().