CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Attributes
EGRegressionModifierV1 Class Reference
Inheritance diagram for EGRegressionModifierV1:
ModifyObjectValueBase

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_
 
CaloGeometry const * caloGeom_
 
edm::ESGetToken< CaloGeometry,
CaloGeometryRecord
caloGeomToken_
 
edm::ESGetToken< GBRForest,
GBRWrapperRcd
condNamesWeight25nsToken_
 
edm::ESGetToken< GBRForest,
GBRWrapperRcd
condNamesWeight50nsToken_
 
EGRegressionModifierCondTokens eleCond25nsTokens_
 
EGRegressionModifierCondTokens eleCond50nsTokens_
 
std::vector< const GBRForestD * > eleForestsMean_
 
std::vector< const GBRForestD * > eleForestsSigma_
 
const GBRForestepForest_
 
int nVtx_
 
EGRegressionModifierCondTokens phoCond25nsTokens_
 
EGRegressionModifierCondTokens phoCond50nsTokens_
 
std::vector< const GBRForestD * > phoForestsMean_
 
std::vector< const GBRForestD * > phoForestsSigma_
 
edm::EDGetTokenT< double > rhoToken_
 
float rhoValue_
 
edm::EDGetTokenT
< reco::VertexCollection
vtxToken_
 

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 61 of file EGRegressionModifierV1.cc.

References edm::ParameterSet::getParameterSet().

62  : ModifyObjectValueBase(conf),
63  eleCond50nsTokens_{conf.getParameterSet("electron_config"), "regressionKey_50ns", "uncertaintyKey_50ns", cc},
64  phoCond50nsTokens_{conf.getParameterSet("photon_config"), "regressionKey_50ns", "uncertaintyKey_50ns", cc},
65  eleCond25nsTokens_{conf.getParameterSet("electron_config"), "regressionKey_25ns", "uncertaintyKey_25ns", cc},
66  phoCond25nsTokens_{conf.getParameterSet("photon_config"), "regressionKey_25ns", "uncertaintyKey_25ns", cc},
67  autoDetectBunchSpacing_(conf.getParameter<bool>("autoDetectBunchSpacing")),
68  bunchspacing_(autoDetectBunchSpacing_ ? 450 : conf.getParameter<int>("manualBunchSpacing")),
69  rhoToken_(cc.consumes(conf.getParameter<edm::InputTag>("rhoCollection"))),
70  vtxToken_(cc.consumes(conf.getParameter<edm::InputTag>("vertexCollection"))),
72  applyExtraHighEnergyProtection_(conf.getParameter<bool>("applyExtraHighEnergyProtection")) {
74  bunchSpacingToken_ = cc.consumes(conf.getParameter<edm::InputTag>("bunchSpacingTag"));
75 
76  auto const& electrons = conf.getParameterSet("electron_config");
77  condNamesWeight50nsToken_ = cc.esConsumes(electrons.getParameter<edm::ESInputTag>("combinationKey_50ns"));
78  condNamesWeight25nsToken_ = cc.esConsumes(electrons.getParameter<edm::ESInputTag>("combinationKey_25ns"));
79 }
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::ESGetToken< GBRForest, GBRWrapperRcd > condNamesWeight25nsToken_
edm::ESGetToken< GBRForest, GBRWrapperRcd > condNamesWeight50nsToken_
EGRegressionModifierCondTokens phoCond50nsTokens_
EGRegressionModifierCondTokens phoCond25nsTokens_
EGRegressionModifierCondTokens eleCond50nsTokens_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
edm::EDGetTokenT< double > rhoToken_
ModifyObjectValueBase(const edm::ParameterSet &conf)
ParameterSet const & getParameterSet(std::string const &) const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
EGRegressionModifierCondTokens eleCond25nsTokens_
edm::EDGetTokenT< unsigned int > bunchSpacingToken_

Member Function Documentation

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

Reimplemented from ModifyObjectValueBase.

Definition at line 103 of file EGRegressionModifierV1.cc.

References funct::abs(), applyExtraHighEnergyProtection_, caloGeom_, reco::GsfElectron::classification(), reco::GsfElectron::correctMomentum(), reco::deltaPhi(), reco::deltaR(), reco::GsfElectron::ecalDriven(), eleForestsMean_, eleForestsSigma_, epForest_, GBRForest::GetResponse(), edm::Ref< C, T, F >::isAvailable(), reco::GsfElectron::isEB(), egammaTools::localEcalClusterCoordsEB(), egammaTools::localEcalClusterCoordsEE(), SiStripPI::max, SiStripPI::mean, nVtx_, reco::GsfElectron::p4(), reco::GsfElectron::setCorrectedEcalEnergy(), reco::GsfElectron::setCorrectedEcalEnergyError(), reco::GsfElectron::showerShape(), mathSSE::sqrt(), reco::GsfElectron::superCluster(), reco::GsfElectron::trackerDrivenSeed(), reco::GsfElectron::trackMomentumAtVtx(), reco::GsfElectron::trackMomentumError(), and histoStyle::weight.

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

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

284  {
285  // regression calculation needs no additional valuemaps
286 
287  std::array<float, 35> eval;
288  const reco::SuperClusterRef& superClus = pho.superCluster();
289  const edm::Ptr<reco::CaloCluster>& theseed = superClus->seed();
290 
291  const int numberOfClusters = superClus->clusters().size();
292  const bool missing_clusters = !superClus->clusters()[numberOfClusters - 1].isAvailable();
293 
294  if (missing_clusters)
295  return; // do not apply corrections in case of missing info (slimmed MiniAOD electrons)
296 
297  const double rawEnergy = superClus->rawEnergy();
298  const auto& ess = pho.showerShapeVariables();
299 
300  // SET INPUTS
301  eval[0] = rawEnergy;
302  eval[1] = pho.r9();
303  eval[2] = superClus->etaWidth();
304  eval[3] = superClus->phiWidth();
305  eval[4] = std::max(0, numberOfClusters - 1);
306  eval[5] = pho.hadronicOverEm();
307  eval[6] = rhoValue_;
308  eval[7] = nVtx_;
309  eval[8] = theseed->eta() - superClus->position().Eta();
310  eval[9] = reco::deltaPhi(theseed->phi(), superClus->position().Phi());
311  eval[10] = theseed->energy() / rawEnergy;
312  eval[11] = ess.e3x3 / ess.e5x5;
313  eval[12] = ess.sigmaIetaIeta;
314  eval[13] = ess.sigmaIphiIphi;
315  eval[14] = ess.sigmaIetaIphi / (ess.sigmaIphiIphi * ess.sigmaIetaIeta);
316  eval[15] = ess.maxEnergyXtal / ess.e5x5;
317  eval[16] = ess.e2nd / ess.e5x5;
318  eval[17] = ess.eTop / ess.e5x5;
319  eval[18] = ess.eBottom / ess.e5x5;
320  eval[19] = ess.eLeft / ess.e5x5;
321  eval[20] = ess.eRight / ess.e5x5;
322  eval[21] = ess.e2x5Max / ess.e5x5;
323  eval[22] = ess.e2x5Left / ess.e5x5;
324  eval[23] = ess.e2x5Right / ess.e5x5;
325  eval[24] = ess.e2x5Top / ess.e5x5;
326  eval[25] = ess.e2x5Bottom / ess.e5x5;
327 
328  const bool isEB = pho.isEB();
329  if (isEB) {
330  EBDetId ebseedid(theseed->seed());
331  eval[26] = pho.e5x5() / theseed->energy();
332  int ieta = ebseedid.ieta();
333  int iphi = ebseedid.iphi();
334  eval[27] = ieta;
335  eval[28] = iphi;
336  int signieta = ieta > 0 ? +1 : -1;
337  eval[29] = (ieta - signieta) % 5;
338  eval[30] = (iphi - 1) % 2;
339  eval[31] = (abs(ieta) <= 25) * ((ieta - signieta)) + (abs(ieta) > 25) * ((ieta - 26 * signieta) % 20);
340  eval[32] = (iphi - 1) % 20;
341  eval[33] = ieta;
342  eval[34] = iphi;
343  } else {
344  EEDetId eeseedid(theseed->seed());
345  eval[26] = superClus->preshowerEnergy() / rawEnergy;
346  eval[27] = superClus->preshowerEnergyPlane1() / rawEnergy;
347  eval[28] = superClus->preshowerEnergyPlane2() / rawEnergy;
348  eval[29] = eeseedid.ix();
349  eval[30] = eeseedid.iy();
350  }
351 
352  //magic numbers for MINUIT-like transformation of BDT output onto limited range
353  //(These should be stored inside the conditions object in the future as well)
354  const double meanlimlow = 0.2;
355  const double meanlimhigh = 2.0;
356  const double meanoffset = meanlimlow + 0.5 * (meanlimhigh - meanlimlow);
357  const double meanscale = 0.5 * (meanlimhigh - meanlimlow);
358 
359  const double sigmalimlow = 0.0002;
360  const double sigmalimhigh = 0.5;
361  const double sigmaoffset = sigmalimlow + 0.5 * (sigmalimhigh - sigmalimlow);
362  const double sigmascale = 0.5 * (sigmalimhigh - sigmalimlow);
363 
364  const int coridx = isEB ? 0 : 1;
365 
366  //these are the actual BDT responses
367  const double rawmean = phoForestsMean_[coridx]->GetResponse(eval.data());
368  const double rawsigma = phoForestsSigma_[coridx]->GetResponse(eval.data());
369  //apply transformation to limited output range (matching the training)
370  const double mean = meanoffset + meanscale * vdt::fast_sin(rawmean);
371  const double sigma = sigmaoffset + sigmascale * vdt::fast_sin(rawsigma);
372 
373  //regression target is ln(Etrue/Eraw)
374  //so corrected energy is ecor=exp(mean)*e, uncertainty is exp(mean)*eraw*sigma=ecor*sigma
375  const double ecor = isEB ? mean * eval[0] : mean * (eval[0] + superClus->preshowerEnergy());
376 
377  const double sigmacor = sigma * ecor;
378  pho.setCorrectedEnergy(reco::Photon::P4type::regression2, ecor, sigmacor, true);
379 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
bool isAvailable() const
Definition: Ref.h:537
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
float e5x5() const
Definition: Photon.h:267
std::vector< const GBRForestD * > phoForestsMean_
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float hadronicOverEm(int depth=0) const
Definition: Photon.h:233
std::vector< const GBRForestD * > phoForestsSigma_
bool isEB() const
Definition: Photon.h:120
float r9() const
Definition: Photon.h:273
const ShowerShape & showerShapeVariables() const
Definition: Photon.h:209
void EGRegressionModifierV1::modifyObject ( pat::Electron ele) const
inlinefinalvirtual

Reimplemented from ModifyObjectValueBase.

Definition at line 29 of file EGRegressionModifierV1.cc.

References modifyObject().

Referenced by modifyObject().

29 { 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 30 of file EGRegressionModifierV1.cc.

References modifyObject().

Referenced by modifyObject().

30 { 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 81 of file EGRegressionModifierV1.cc.

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

81  {
84  }
85  rhoValue_ = evt.get(rhoToken_);
86  nVtx_ = evt.get(vtxToken_).size();
87 }
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
edm::EDGetTokenT< double > rhoToken_
edm::EDGetTokenT< unsigned int > bunchSpacingToken_
void EGRegressionModifierV1::setEventContent ( const edm::EventSetup evs)
finalvirtual

Reimplemented from ModifyObjectValueBase.

Definition at line 89 of file EGRegressionModifierV1.cc.

References bunchspacing_, caloGeom_, caloGeomToken_, condNamesWeight25nsToken_, condNamesWeight50nsToken_, eleCond25nsTokens_, eleCond50nsTokens_, eleForestsMean_, eleForestsSigma_, epForest_, edm::EventSetup::getData(), EGRegressionModifierCondTokens::mean, phoCond25nsTokens_, phoCond50nsTokens_, phoForestsMean_, phoForestsSigma_, retrieveGBRForests(), and EGRegressionModifierCondTokens::sigma.

89  {
91 
95 
99 
101 }
edm::ESGetToken< GBRForest, GBRWrapperRcd > condNamesWeight25nsToken_
edm::ESGetToken< GBRForest, GBRWrapperRcd > condNamesWeight50nsToken_
EGRegressionModifierCondTokens phoCond50nsTokens_
std::vector< edm::ESGetToken< GBRForestD, GBRDWrapperRcd > > mean
EGRegressionModifierCondTokens phoCond25nsTokens_
EGRegressionModifierCondTokens eleCond50nsTokens_
CaloGeometry const * caloGeom_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:128
std::vector< const GBRForestD * > phoForestsMean_
std::vector< const GBRForestD * > eleForestsSigma_
std::vector< const GBRForestD * > eleForestsMean_
std::vector< const GBRForestD * > phoForestsSigma_
std::vector< edm::ESGetToken< GBRForestD, GBRDWrapperRcd > > sigma
std::vector< const GBRForestD * > retrieveGBRForests(edm::EventSetup const &evs, std::vector< edm::ESGetToken< GBRForestD, GBRDWrapperRcd >> const &tokens)
EGRegressionModifierCondTokens eleCond25nsTokens_

Member Data Documentation

const bool EGRegressionModifierV1::applyExtraHighEnergyProtection_
private

Definition at line 50 of file EGRegressionModifierV1.cc.

Referenced by modifyObject().

const bool EGRegressionModifierV1::autoDetectBunchSpacing_
private

Definition at line 41 of file EGRegressionModifierV1.cc.

Referenced by setEvent().

int EGRegressionModifierV1::bunchspacing_
private

Definition at line 42 of file EGRegressionModifierV1.cc.

Referenced by setEvent(), and setEventContent().

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

Definition at line 43 of file EGRegressionModifierV1.cc.

Referenced by setEvent().

CaloGeometry const* EGRegressionModifierV1::caloGeom_
private

Definition at line 49 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

edm::ESGetToken<CaloGeometry, CaloGeometryRecord> EGRegressionModifierV1::caloGeomToken_
private

Definition at line 48 of file EGRegressionModifierV1.cc.

Referenced by setEventContent().

edm::ESGetToken<GBRForest, GBRWrapperRcd> EGRegressionModifierV1::condNamesWeight25nsToken_
private

Definition at line 39 of file EGRegressionModifierV1.cc.

Referenced by setEventContent().

edm::ESGetToken<GBRForest, GBRWrapperRcd> EGRegressionModifierV1::condNamesWeight50nsToken_
private

Definition at line 38 of file EGRegressionModifierV1.cc.

Referenced by setEventContent().

EGRegressionModifierCondTokens EGRegressionModifierV1::eleCond25nsTokens_
private

Definition at line 35 of file EGRegressionModifierV1.cc.

Referenced by setEventContent().

EGRegressionModifierCondTokens EGRegressionModifierV1::eleCond50nsTokens_
private

Definition at line 33 of file EGRegressionModifierV1.cc.

Referenced by setEventContent().

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

Definition at line 54 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

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

Definition at line 55 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

const GBRForest* EGRegressionModifierV1::epForest_
private

Definition at line 56 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

int EGRegressionModifierV1::nVtx_
private

Definition at line 46 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEvent().

EGRegressionModifierCondTokens EGRegressionModifierV1::phoCond25nsTokens_
private

Definition at line 36 of file EGRegressionModifierV1.cc.

Referenced by setEventContent().

EGRegressionModifierCondTokens EGRegressionModifierV1::phoCond50nsTokens_
private

Definition at line 34 of file EGRegressionModifierV1.cc.

Referenced by setEventContent().

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

Definition at line 52 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

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

Definition at line 53 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEventContent().

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

Definition at line 45 of file EGRegressionModifierV1.cc.

Referenced by setEvent().

float EGRegressionModifierV1::rhoValue_
private

Definition at line 44 of file EGRegressionModifierV1.cc.

Referenced by modifyObject(), and setEvent().

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

Definition at line 47 of file EGRegressionModifierV1.cc.

Referenced by setEvent().