CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes
EGEnergyCorrector Class Reference

#include <EGEnergyCorrector.h>

Public Member Functions

std::pair< double, double > CorrectedEnergyWithError (const reco::Photon &p, const reco::VertexCollection &vtxcol, EcalClusterLazyTools &clustertools, CaloGeometry const &caloGeometry)
 
std::pair< double, double > CorrectedEnergyWithErrorV3 (const reco::Photon &p, const reco::VertexCollection &vtxcol, double rho, EcalClusterLazyTools &clustertools, CaloGeometry const &caloGeometry, bool applyRescale=false)
 
void Initialize (const edm::EventSetup &iSetup, std::string regweights, bool weightsFromDB=false)
 
bool IsInitialized () const
 
 ~EGEnergyCorrector ()
 

Protected Attributes

bool fIsInitialized = false
 
bool fOwnsForests = false
 
const GBRForestfReadereb = nullptr
 
const GBRForestfReaderebvariance = nullptr
 
const GBRForestfReaderee = nullptr
 
const GBRForestfReadereevariance = nullptr
 
std::array< float, 73 > fVals
 

Detailed Description

Definition at line 22 of file EGEnergyCorrector.h.

Constructor & Destructor Documentation

◆ ~EGEnergyCorrector()

EGEnergyCorrector::~EGEnergyCorrector ( )

Definition at line 17 of file EGEnergyCorrector.cc.

17  {
18  if (fOwnsForests) {
19  if (fReadereb)
20  delete fReadereb;
22  delete fReaderebvariance;
23  if (fReaderee)
24  delete fReaderee;
26  delete fReadereevariance;
27  }
28 }

Member Function Documentation

◆ CorrectedEnergyWithError()

std::pair< double, double > EGEnergyCorrector::CorrectedEnergyWithError ( const reco::Photon p,
const reco::VertexCollection vtxcol,
EcalClusterLazyTools clustertools,
CaloGeometry const &  caloGeometry 
)

Definition at line 77 of file EGEnergyCorrector.cc.

80  {
81  const SuperClusterRef s = p.superCluster();
82  const CaloClusterPtr b = s->seed(); //seed basic cluster
83 
84  //highest energy basic cluster excluding seed basic cluster
86  Double_t ebcmax = -99.;
87  for (reco::CaloCluster_iterator bit = s->clustersBegin(); bit != s->clustersEnd(); ++bit) {
88  const CaloClusterPtr bc = *bit;
89  if (bc->energy() > ebcmax && bc != b) {
90  b2 = bc;
91  ebcmax = bc->energy();
92  }
93  }
94 
95  //lowest energy basic cluster excluding seed (for pileup mitigation)
96  CaloClusterPtr bclast;
97  Double_t ebcmin = 1e6;
98  for (reco::CaloCluster_iterator bit = s->clustersBegin(); bit != s->clustersEnd(); ++bit) {
99  const CaloClusterPtr bc = *bit;
100  if (bc->energy() < ebcmin && bc != b) {
101  bclast = bc;
102  ebcmin = bc->energy();
103  }
104  }
105 
106  //2nd lowest energy basic cluster excluding seed (for pileup mitigation)
107  CaloClusterPtr bclast2;
108  ebcmin = 1e6;
109  for (reco::CaloCluster_iterator bit = s->clustersBegin(); bit != s->clustersEnd(); ++bit) {
110  const CaloClusterPtr bc = *bit;
111  if (bc->energy() < ebcmin && bc != b && bc != bclast) {
112  bclast2 = bc;
113  ebcmin = bc->energy();
114  }
115  }
116 
117  Bool_t isbarrel = b->hitsAndFractions().at(0).first.subdetId() == EcalBarrel;
118  Bool_t hasbc2 = b2.isNonnull() && b2->energy() > 0.;
119  Bool_t hasbclast = bclast.isNonnull() && bclast->energy() > 0.;
120  Bool_t hasbclast2 = bclast2.isNonnull() && bclast2->energy() > 0.;
121 
122  if (isbarrel) {
123  //basic supercluster variables
124  fVals[0] = s->rawEnergy();
125  fVals[1] = p.r9();
126  fVals[2] = s->eta();
127  fVals[3] = s->phi();
128  fVals[4] = p.e5x5() / s->rawEnergy();
129  fVals[5] = p.hadronicOverEm();
130  fVals[6] = s->etaWidth();
131  fVals[7] = s->phiWidth();
132 
133  //seed basic cluster variables
134  double bemax = clustertools.eMax(*b);
135  double be2nd = clustertools.e2nd(*b);
136  double betop = clustertools.eTop(*b);
137  double bebottom = clustertools.eBottom(*b);
138  double beleft = clustertools.eLeft(*b);
139  double beright = clustertools.eRight(*b);
140 
141  fVals[8] = b->eta() - s->eta();
142  fVals[9] = reco::deltaPhi(b->phi(), s->phi());
143  fVals[10] = b->energy() / s->rawEnergy();
144  fVals[11] = clustertools.e3x3(*b) / b->energy();
145  fVals[12] = clustertools.e5x5(*b) / b->energy();
146  fVals[13] = sqrt(clustertools.localCovariances(*b)[0]); //sigietaieta
147  fVals[14] = sqrt(clustertools.localCovariances(*b)[2]); //sigiphiiphi
148  fVals[15] = clustertools.localCovariances(*b)[1]; //sigietaiphi
149  fVals[16] = bemax / b->energy(); //crystal energy ratio gap variables
150  fVals[17] = log(be2nd / bemax);
151  fVals[18] = log(betop / bemax);
152  fVals[19] = log(bebottom / bemax);
153  fVals[20] = log(beleft / bemax);
154  fVals[21] = log(beright / bemax);
155  fVals[22] = (betop - bebottom) / (betop + bebottom);
156  fVals[23] = (beleft - beright) / (beleft + beright);
157 
158  double bc2emax = hasbc2 ? clustertools.eMax(*b2) : 0.;
159  double bc2e2nd = hasbc2 ? clustertools.e2nd(*b2) : 0.;
160  double bc2etop = hasbc2 ? clustertools.eTop(*b2) : 0.;
161  double bc2ebottom = hasbc2 ? clustertools.eBottom(*b2) : 0.;
162  double bc2eleft = hasbc2 ? clustertools.eLeft(*b2) : 0.;
163  double bc2eright = hasbc2 ? clustertools.eRight(*b2) : 0.;
164 
165  fVals[24] = hasbc2 ? (b2->eta() - s->eta()) : 0.;
166  fVals[25] = hasbc2 ? reco::deltaPhi(b2->phi(), s->phi()) : 0.;
167  fVals[26] = hasbc2 ? b2->energy() / s->rawEnergy() : 0.;
168  fVals[27] = hasbc2 ? clustertools.e3x3(*b2) / b2->energy() : 0.;
169  fVals[28] = hasbc2 ? clustertools.e5x5(*b2) / b2->energy() : 0.;
170  fVals[29] = hasbc2 ? sqrt(clustertools.localCovariances(*b2)[0]) : 0.;
171  fVals[30] = hasbc2 ? sqrt(clustertools.localCovariances(*b2)[2]) : 0.;
172  fVals[31] = hasbc2 ? clustertools.localCovariances(*b)[1] : 0.;
173  fVals[32] = hasbc2 ? bc2emax / b2->energy() : 0.;
174  fVals[33] = hasbc2 ? log(bc2e2nd / bc2emax) : 0.;
175  fVals[34] = hasbc2 ? log(bc2etop / bc2emax) : 0.;
176  fVals[35] = hasbc2 ? log(bc2ebottom / bc2emax) : 0.;
177  fVals[36] = hasbc2 ? log(bc2eleft / bc2emax) : 0.;
178  fVals[37] = hasbc2 ? log(bc2eright / bc2emax) : 0.;
179  fVals[38] = hasbc2 ? (bc2etop - bc2ebottom) / (bc2etop + bc2ebottom) : 0.;
180  fVals[39] = hasbc2 ? (bc2eleft - bc2eright) / (bc2eleft + bc2eright) : 0.;
181 
182  fVals[40] = hasbclast ? (bclast->eta() - s->eta()) : 0.;
183  fVals[41] = hasbclast ? reco::deltaPhi(bclast->phi(), s->phi()) : 0.;
184  fVals[42] = hasbclast ? bclast->energy() / s->rawEnergy() : 0.;
185  fVals[43] = hasbclast ? clustertools.e3x3(*bclast) / bclast->energy() : 0.;
186  fVals[44] = hasbclast ? clustertools.e5x5(*bclast) / bclast->energy() : 0.;
187  fVals[45] = hasbclast ? sqrt(clustertools.localCovariances(*bclast)[0]) : 0.;
188  fVals[46] = hasbclast ? sqrt(clustertools.localCovariances(*bclast)[2]) : 0.;
189  fVals[47] = hasbclast ? clustertools.localCovariances(*bclast)[1] : 0.;
190 
191  fVals[48] = hasbclast2 ? (bclast2->eta() - s->eta()) : 0.;
192  fVals[49] = hasbclast2 ? reco::deltaPhi(bclast2->phi(), s->phi()) : 0.;
193  fVals[50] = hasbclast2 ? bclast2->energy() / s->rawEnergy() : 0.;
194  fVals[51] = hasbclast2 ? clustertools.e3x3(*bclast2) / bclast2->energy() : 0.;
195  fVals[52] = hasbclast2 ? clustertools.e5x5(*bclast2) / bclast2->energy() : 0.;
196  fVals[53] = hasbclast2 ? sqrt(clustertools.localCovariances(*bclast2)[0]) : 0.;
197  fVals[54] = hasbclast2 ? sqrt(clustertools.localCovariances(*bclast2)[2]) : 0.;
198  fVals[55] = hasbclast2 ? clustertools.localCovariances(*bclast2)[1] : 0.;
199 
200  //local coordinates and crystal indices
201 
202  //seed cluster
203  float betacry, bphicry, bthetatilt, bphitilt;
204  int bieta, biphi;
205  egammaTools::localEcalClusterCoordsEB(*b, caloGeometry, betacry, bphicry, bieta, biphi, bthetatilt, bphitilt);
206 
207  fVals[56] = bieta; //crystal ieta
208  fVals[57] = biphi; //crystal iphi
209  fVals[58] = bieta % 5; //submodule boundary eta symmetry
210  fVals[59] = biphi % 2; //submodule boundary phi symmetry
211  fVals[60] = (TMath::Abs(bieta) <= 25) * (bieta % 25) +
212  (TMath::Abs(bieta) > 25) *
213  ((bieta - 25 * TMath::Abs(bieta) / bieta) % 20); //module boundary eta approximate symmetry
214  fVals[61] = biphi % 20; //module boundary phi symmetry
215  fVals[62] = betacry; //local coordinates with respect to closest crystal center at nominal shower depth
216  fVals[63] = bphicry;
217 
218  //2nd cluster (meaningful gap corrections for converted photons)
219  float bc2etacry, bc2phicry, bc2thetatilt, bc2phitilt;
220  int bc2ieta, bc2iphi;
221  if (hasbc2)
223  *b2, caloGeometry, bc2etacry, bc2phicry, bc2ieta, bc2iphi, bc2thetatilt, bc2phitilt);
224 
225  fVals[64] = hasbc2 ? bc2ieta : 0.;
226  fVals[65] = hasbc2 ? bc2iphi : 0.;
227  fVals[66] = hasbc2 ? bc2ieta % 5 : 0.;
228  fVals[67] = hasbc2 ? bc2iphi % 2 : 0.;
229  fVals[68] = hasbc2 ? (TMath::Abs(bc2ieta) <= 25) * (bc2ieta % 25) +
230  (TMath::Abs(bc2ieta) > 25) * ((bc2ieta - 25 * TMath::Abs(bc2ieta) / bc2ieta) % 20)
231  : 0.;
232  fVals[69] = hasbc2 ? bc2iphi % 20 : 0.;
233  fVals[70] = hasbc2 ? bc2etacry : 0.;
234  fVals[71] = hasbc2 ? bc2phicry : 0.;
235 
236  fVals[72] = vtxcol.size();
237 
238  } else {
239  fVals[0] = s->rawEnergy();
240  fVals[1] = p.r9();
241  fVals[2] = s->eta();
242  fVals[3] = s->phi();
243  fVals[4] = p.e5x5() / s->rawEnergy();
244  fVals[5] = s->etaWidth();
245  fVals[6] = s->phiWidth();
246  fVals[7] = vtxcol.size();
247  }
248 
249  const Double_t varscale = 1.253;
250  Double_t den;
251  const GBRForest *reader;
252  const GBRForest *readervar;
253  if (isbarrel) {
254  den = s->rawEnergy();
255  reader = fReadereb;
256  readervar = fReaderebvariance;
257  } else {
258  den = s->rawEnergy() + s->preshowerEnergy();
259  reader = fReaderee;
260  readervar = fReadereevariance;
261  }
262 
263  Double_t ecor = reader->GetResponse(fVals.data()) * den;
264  Double_t ecorerr = readervar->GetResponse(fVals.data()) * den * varscale;
265 
266  //printf("ecor = %5f, ecorerr = %5f\n",ecor,ecorerr);
267 
268  return {ecor, ecorerr};
269 }

References Abs(), b, b2, reco::deltaPhi(), EcalBarrel, reco::CaloCluster::energy(), reco::CaloCluster::eta(), GBRForest::GetResponse(), edm::Ptr< T >::isNonnull(), egammaTools::localEcalClusterCoordsEB(), dqm-mbProfile::log, AlCaHLTBitMon_ParallelJobs::p, reco::CaloCluster::phi(), DQM::reader, alignCSCRings::s, and mathSSE::sqrt().

Referenced by EGEnergyAnalyzer::analyze().

◆ CorrectedEnergyWithErrorV3()

std::pair< double, double > EGEnergyCorrector::CorrectedEnergyWithErrorV3 ( const reco::Photon p,
const reco::VertexCollection vtxcol,
double  rho,
EcalClusterLazyTools clustertools,
CaloGeometry const &  caloGeometry,
bool  applyRescale = false 
)

Definition at line 272 of file EGEnergyCorrector.cc.

277  {
278  const SuperClusterRef s = p.superCluster();
279  const CaloClusterPtr b = s->seed(); //seed basic cluster
280 
281  Bool_t isbarrel = b->hitsAndFractions().at(0).first.subdetId() == EcalBarrel;
282 
283  //basic supercluster variables
284  fVals[0] = s->rawEnergy();
285  fVals[1] = s->eta();
286  fVals[2] = s->phi();
287  fVals[3] = p.r9();
288  fVals[4] = p.e5x5() / s->rawEnergy();
289  fVals[5] = s->etaWidth();
290  fVals[6] = s->phiWidth();
291  fVals[7] = s->clustersSize();
292  fVals[8] = p.hadTowOverEm();
293  fVals[9] = rho;
294  fVals[10] = vtxcol.size();
295 
296  //seed basic cluster variables
297  double bemax = clustertools.eMax(*b);
298  double be2nd = clustertools.e2nd(*b);
299  double betop = clustertools.eTop(*b);
300  double bebottom = clustertools.eBottom(*b);
301  double beleft = clustertools.eLeft(*b);
302  double beright = clustertools.eRight(*b);
303 
304  double be2x5max = clustertools.e2x5Max(*b);
305  double be2x5top = clustertools.e2x5Top(*b);
306  double be2x5bottom = clustertools.e2x5Bottom(*b);
307  double be2x5left = clustertools.e2x5Left(*b);
308  double be2x5right = clustertools.e2x5Right(*b);
309 
310  fVals[11] = b->eta() - s->eta();
311  fVals[12] = reco::deltaPhi(b->phi(), s->phi());
312  fVals[13] = b->energy() / s->rawEnergy();
313  fVals[14] = clustertools.e3x3(*b) / b->energy();
314  fVals[15] = clustertools.e5x5(*b) / b->energy();
315  fVals[16] = sqrt(clustertools.localCovariances(*b)[0]); //sigietaieta
316  fVals[17] = sqrt(clustertools.localCovariances(*b)[2]); //sigiphiiphi
317  fVals[18] = clustertools.localCovariances(*b)[1]; //sigietaiphi
318  fVals[19] = bemax / b->energy(); //crystal energy ratio gap variables
319  fVals[20] = be2nd / b->energy();
320  fVals[21] = betop / b->energy();
321  fVals[22] = bebottom / b->energy();
322  fVals[23] = beleft / b->energy();
323  fVals[24] = beright / b->energy();
324  fVals[25] = be2x5max / b->energy(); //crystal energy ratio gap variables
325  fVals[26] = be2x5top / b->energy();
326  fVals[27] = be2x5bottom / b->energy();
327  fVals[28] = be2x5left / b->energy();
328  fVals[29] = be2x5right / b->energy();
329 
330  if (isbarrel) {
331  //local coordinates and crystal indices (barrel only)
332 
333  //seed cluster
334  float betacry, bphicry, bthetatilt, bphitilt;
335  int bieta, biphi;
336  egammaTools::localEcalClusterCoordsEB(*b, caloGeometry, betacry, bphicry, bieta, biphi, bthetatilt, bphitilt);
337 
338  fVals[30] = bieta; //crystal ieta
339  fVals[31] = biphi; //crystal iphi
340  fVals[32] = bieta % 5; //submodule boundary eta symmetry
341  fVals[33] = biphi % 2; //submodule boundary phi symmetry
342  fVals[34] = (TMath::Abs(bieta) <= 25) * (bieta % 25) +
343  (TMath::Abs(bieta) > 25) *
344  ((bieta - 25 * TMath::Abs(bieta) / bieta) % 20); //module boundary eta approximate symmetry
345  fVals[35] = biphi % 20; //module boundary phi symmetry
346  fVals[36] = betacry; //local coordinates with respect to closest crystal center at nominal shower depth
347  fVals[37] = bphicry;
348 
349  } else {
350  //preshower energy ratio (endcap only)
351  fVals[30] = s->preshowerEnergy() / s->rawEnergy();
352  }
353 
354  // if (isbarrel) {
355  // for (int i=0; i<38; ++i) printf("%i: %5f\n",i,fVals[i]);
356  // }
357  // else for (int i=0; i<31; ++i) printf("%i: %5f\n",i,fVals[i]);
358 
359  Double_t den;
360  const GBRForest *reader;
361  const GBRForest *readervar;
362  if (isbarrel) {
363  den = s->rawEnergy();
364  reader = fReadereb;
365  readervar = fReaderebvariance;
366  } else {
367  den = s->rawEnergy() + s->preshowerEnergy();
368  reader = fReaderee;
369  readervar = fReadereevariance;
370  }
371 
372  Double_t ecor = reader->GetResponse(fVals.data()) * den;
373 
374  //apply shower shape rescaling - for Monte Carlo only, and only for calculation of energy uncertainty
375  if (applyRescale) {
376  if (isbarrel) {
377  fVals[3] = 1.0045 * p.r9() + 0.001; //r9
378  fVals[5] = 1.04302 * s->etaWidth() - 0.000618; //etawidth
379  fVals[6] = 1.00002 * s->phiWidth() - 0.000371; //phiwidth
380  fVals[14] = fVals[3] * s->rawEnergy() / b->energy(); //compute consistent e3x3/eseed after r9 rescaling
381  if (fVals[15] <= 1.0) // rescale e5x5/eseed only if value is <=1.0, don't allow scaled values to exceed 1.0
382  fVals[15] = TMath::Min(1.0, 1.0022 * p.e5x5() / b->energy());
383 
384  fVals[4] =
385  fVals[15] * b->energy() / s->rawEnergy(); // compute consistent e5x5()/rawEnergy() after e5x5/eseed resacling
386 
387  fVals[16] = 0.891832 * sqrt(clustertools.localCovariances(*b)[0]) + 0.0009133; //sigietaieta
388  fVals[17] = 0.993 * sqrt(clustertools.localCovariances(*b)[2]); //sigiphiiphi
389 
390  fVals[19] = 1.012 * bemax / b->energy(); //crystal energy ratio gap variables
391  fVals[20] = 1.0 * be2nd / b->energy();
392  fVals[21] = 0.94 * betop / b->energy();
393  fVals[22] = 0.94 * bebottom / b->energy();
394  fVals[23] = 0.94 * beleft / b->energy();
395  fVals[24] = 0.94 * beright / b->energy();
396  fVals[25] = 1.006 * be2x5max / b->energy(); //crystal energy ratio gap variables
397  fVals[26] = 1.09 * be2x5top / b->energy();
398  fVals[27] = 1.09 * be2x5bottom / b->energy();
399  fVals[28] = 1.09 * be2x5left / b->energy();
400  fVals[29] = 1.09 * be2x5right / b->energy();
401 
402  } else {
403  fVals[3] = 1.0086 * p.r9() - 0.0007; //r9
404  fVals[4] = TMath::Min(1.0, 1.0022 * p.e5x5() / s->rawEnergy()); //e5x5/rawenergy
405  fVals[5] = 0.903254 * s->etaWidth() + 0.001346; //etawidth
406  fVals[6] = 0.99992 * s->phiWidth() + 4.8e-07; //phiwidth
407  fVals[13] =
408  TMath::Min(1.0, 1.0022 * b->energy() / s->rawEnergy()); //eseed/rawenergy (practically equivalent to e5x5)
409 
410  fVals[14] = fVals[3] * s->rawEnergy() / b->energy(); //compute consistent e3x3/eseed after r9 rescaling
411 
412  fVals[16] = 0.9947 * sqrt(clustertools.localCovariances(*b)[0]) + 0.00003; //sigietaieta
413 
414  fVals[19] = 1.005 * bemax / b->energy(); //crystal energy ratio gap variables
415  fVals[20] = 1.02 * be2nd / b->energy();
416  fVals[21] = 0.96 * betop / b->energy();
417  fVals[22] = 0.96 * bebottom / b->energy();
418  fVals[23] = 0.96 * beleft / b->energy();
419  fVals[24] = 0.96 * beright / b->energy();
420  fVals[25] = 1.0075 * be2x5max / b->energy(); //crystal energy ratio gap variables
421  fVals[26] = 1.13 * be2x5top / b->energy();
422  fVals[27] = 1.13 * be2x5bottom / b->energy();
423  fVals[28] = 1.13 * be2x5left / b->energy();
424  fVals[29] = 1.13 * be2x5right / b->energy();
425  }
426  }
427 
428  Double_t ecorerr = readervar->GetResponse(fVals.data()) * den;
429 
430  //printf("ecor = %5f, ecorerr = %5f\n",ecor,ecorerr);
431 
432  return {ecor, ecorerr};
433 }

References Abs(), b, reco::deltaPhi(), EcalBarrel, GBRForest::GetResponse(), egammaTools::localEcalClusterCoordsEB(), Min(), AlCaHLTBitMon_ParallelJobs::p, DQM::reader, alignCSCRings::s, and mathSSE::sqrt().

◆ Initialize()

void EGEnergyCorrector::Initialize ( const edm::EventSetup iSetup,
std::string  regweights,
bool  weightsFromDB = false 
)

Definition at line 31 of file EGEnergyCorrector.cc.

31  {
32  fIsInitialized = true;
33 
34  if (fOwnsForests) {
35  if (fReadereb)
36  delete fReadereb;
38  delete fReaderebvariance;
39  if (fReaderee)
40  delete fReaderee;
42  delete fReadereevariance;
43  }
44 
45  fVals.fill(0.0f);
46 
47  if (weightsFromDB) { //weights from event setup
48 
49  edm::ESHandle<GBRForest> readereb;
50  edm::ESHandle<GBRForest> readerebvar;
51  edm::ESHandle<GBRForest> readeree;
52  edm::ESHandle<GBRForest> readereevar;
53 
54  iSetup.get<GBRWrapperRcd>().get(regweights + "_EBCorrection", readereb);
55  iSetup.get<GBRWrapperRcd>().get(regweights + "_EBUncertainty", readerebvar);
56  iSetup.get<GBRWrapperRcd>().get(regweights + "_EECorrection", readeree);
57  iSetup.get<GBRWrapperRcd>().get(regweights + "_EEUncertainty", readereevar);
58 
59  fReadereb = readereb.product();
60  fReaderebvariance = readerebvar.product();
61  fReaderee = readeree.product();
62  fReadereevariance = readereevar.product();
63 
64  } else { //weights from root file
65  fOwnsForests = true;
66 
67  TFile *fgbr = TFile::Open(regweights.c_str(), "READ");
68  fReadereb = (GBRForest *)fgbr->Get("EBCorrection");
69  fReaderebvariance = (GBRForest *)fgbr->Get("EBUncertainty");
70  fReaderee = (GBRForest *)fgbr->Get("EECorrection");
71  fReadereevariance = (GBRForest *)fgbr->Get("EEUncertainty");
72  fgbr->Close();
73  }
74 }

References f, reco::get(), edm::EventSetup::get(), and edm::ESHandle< T >::product().

Referenced by EGEnergyAnalyzer::analyze().

◆ IsInitialized()

bool EGEnergyCorrector::IsInitialized ( ) const
inline

Definition at line 27 of file EGEnergyCorrector.h.

27 { return fIsInitialized; }

References fIsInitialized.

Referenced by EGEnergyAnalyzer::analyze().

Member Data Documentation

◆ fIsInitialized

bool EGEnergyCorrector::fIsInitialized = false
protected

Definition at line 47 of file EGEnergyCorrector.h.

Referenced by IsInitialized().

◆ fOwnsForests

bool EGEnergyCorrector::fOwnsForests = false
protected

Definition at line 48 of file EGEnergyCorrector.h.

◆ fReadereb

const GBRForest* EGEnergyCorrector::fReadereb = nullptr
protected

Definition at line 42 of file EGEnergyCorrector.h.

◆ fReaderebvariance

const GBRForest* EGEnergyCorrector::fReaderebvariance = nullptr
protected

Definition at line 43 of file EGEnergyCorrector.h.

◆ fReaderee

const GBRForest* EGEnergyCorrector::fReaderee = nullptr
protected

Definition at line 44 of file EGEnergyCorrector.h.

◆ fReadereevariance

const GBRForest* EGEnergyCorrector::fReadereevariance = nullptr
protected

Definition at line 45 of file EGEnergyCorrector.h.

◆ fVals

std::array<float, 73> EGEnergyCorrector::fVals
protected

Definition at line 49 of file EGEnergyCorrector.h.

reco::CaloCluster::phi
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
EGEnergyCorrector::fReadereb
const GBRForest * fReadereb
Definition: EGEnergyCorrector.h:42
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
reco::deltaPhi
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
GBRForest
Definition: GBRForest.h:25
EGEnergyCorrector::fOwnsForests
bool fOwnsForests
Definition: EGEnergyCorrector.h:48
edm::PtrVectorItr
Definition: PtrVector.h:51
EGEnergyCorrector::fReaderebvariance
const GBRForest * fReaderebvariance
Definition: EGEnergyCorrector.h:43
b2
static constexpr float b2
Definition: L1EGammaCrystalsEmulatorProducer.cc:83
EcalBarrel
Definition: EcalSubdetector.h:10
edm::Ref< SuperClusterCollection >
GBRForest::GetResponse
double GetResponse(const float *vector) const
Definition: GBRForest.h:49
alignCSCRings.s
s
Definition: alignCSCRings.py:92
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
Abs
T Abs(T a)
Definition: MathUtil.h:49
DQM.reader
reader
Definition: DQM.py:105
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
edm::ESHandle< GBRForest >
b
double b
Definition: hdecay.h:118
DDAxes::rho
egammaTools::localEcalClusterCoordsEB
void localEcalClusterCoordsEB(const reco::CaloCluster &bclus, const CaloGeometry &geom, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt)
Definition: EcalClusterLocal.cc:14
reco::CaloCluster::eta
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
EGEnergyCorrector::fIsInitialized
bool fIsInitialized
Definition: EGEnergyCorrector.h:47
EGEnergyCorrector::fVals
std::array< float, 73 > fVals
Definition: EGEnergyCorrector.h:49
get
#define get
edm::Ptr< CaloCluster >
GBRWrapperRcd
Definition: GBRWrapperRcd.h:24
edm::Ptr::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
EGEnergyCorrector::fReadereevariance
const GBRForest * fReadereevariance
Definition: EGEnergyCorrector.h:45
Min
T Min(T a, T b)
Definition: MathUtil.h:39
reco::CaloCluster::energy
double energy() const
cluster energy
Definition: CaloCluster.h:149
EGEnergyCorrector::fReaderee
const GBRForest * fReaderee
Definition: EGEnergyCorrector.h:44