CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
ElectronEnergyCorrector Class Reference

#include <ElectronEnergyCorrector.h>

Public Member Functions

void classBasedParameterizationEnergy (reco::GsfElectron &, const reco::BeamSpot &bs)
 
void classBasedParameterizationUncertainty (reco::GsfElectron &)
 
 ElectronEnergyCorrector (EcalClusterFunctionBaseClass *crackCorrectionFunction)
 
void simpleParameterizationUncertainty (reco::GsfElectron &)
 

Private Member Functions

float fBremEta (float sigmaPhiSigmaEta, float eta, int algorithm, reco::GsfElectron::Classification cl) const
 
float fEnergy (float e, int algorithm, reco::GsfElectron::Classification cl) const
 
float fEt (float et, int algorithm, reco::GsfElectron::Classification cl) const
 
float fEta (float energy, float eta, int algorithm) const
 
double fEtaBarrelBad (double scEta) const
 
double fEtaBarrelGood (double scEta) const
 
double fEtaEndcapBad (double scEta) const
 
double fEtaEndcapGood (double scEta) const
 

Private Attributes

EcalClusterFunctionBaseClasscrackCorrectionFunction_
 

Detailed Description

Definition at line 10 of file ElectronEnergyCorrector.h.

Constructor & Destructor Documentation

ElectronEnergyCorrector::ElectronEnergyCorrector ( EcalClusterFunctionBaseClass crackCorrectionFunction)
inline

Definition at line 14 of file ElectronEnergyCorrector.h.

15  : crackCorrectionFunction_(crackCorrectionFunction) {}
EcalClusterFunctionBaseClass * crackCorrectionFunction_

Member Function Documentation

void ElectronEnergyCorrector::classBasedParameterizationEnergy ( reco::GsfElectron electron,
const reco::BeamSpot bs 
)

Definition at line 56 of file ElectronEnergyCorrector.cc.

References reco::GsfElectron::classification(), alignCSCRings::corr, relval_parameters_module::energy, reco::GsfElectron::GAP, reco::GsfElectron::isEB(), reco::GsfElectron::isEcalEnergyCorrected(), reco::GsfElectron::isEE(), reco::return(), reco::GsfElectron::setCorrectedEcalEnergy(), reco::GsfElectron::superCluster(), and reco::GsfElectron::UNKNOWN.

Referenced by GsfElectronAlgo::createElectron().

57  {
58  if (electron.isEcalEnergyCorrected())
59  {
60  edm::LogWarning("ElectronEnergyCorrector::classBasedParameterizationEnergy")<<"already done" ;
61  return ;
62  }
63 
65  if ( (elClass <= reco::GsfElectron::UNKNOWN) ||
66  (elClass>reco::GsfElectron::GAP) )
67  {
68  edm::LogWarning("ElectronEnergyCorrector::classBasedParameterizationEnergy")<<"unexpected classification" ;
69  return ;
70  }
71 
72  // new corrections from N. Chanon et al., taken from EcalClusterCorrectionObjectSpecific.cc
73  float corr = 1.;
74  float corr2 = 1.;
75  float energy = electron.superCluster()->energy() ;
76  float newEnergy = energy;
77 
78  //int subdet = electron.superCluster()->seed()->hitsAndFractions()[0].first.subdetId();
79 
80  if (electron.isEB())
81  {
82  float cetacorr = fEta(electron.superCluster()->rawEnergy(), electron.superCluster()->eta(), 0)/electron.superCluster()->rawEnergy();
83  energy = electron.superCluster()->rawEnergy()*cetacorr; //previously in CMSSW
84  //energy = superCluster.rawEnergy()*fEta(e5x5, superCluster.seed()->eta(), 0)/e5x5;
85  }
86  else if (electron.isEE())
87  {
88  energy = electron.superCluster()->rawEnergy()+electron.superCluster()->preshowerEnergy();
89  }
90  else
91  { edm::LogWarning("ElectronEnergyCorrector::classBasedParameterizationEnergy")<<"nor barrel neither endcap electron !" ; }
92 
93  corr = fBremEta(electron.superCluster()->phiWidth()/electron.superCluster()->etaWidth(), electron.superCluster()->eta(), 0,elClass);
94 
95  float et = energy*TMath::Sin(2*TMath::ATan(TMath::Exp(-electron.superCluster()->eta())))/corr;
96 
97  if (electron.isEB()) { corr2 = corr * fEt(et, 0,elClass) ; }
98  else if (electron.isEE()) { corr2 = corr * fEnergy(energy/corr, 1,elClass) ; }
99  else { edm::LogWarning("ElectronEnergyCorrector::classBasedParameterizationEnergy")<<"nor barrel neither endcap electron !" ; }
100 
101  newEnergy = energy/corr2;
102 
103  // cracks
104  double crackcor = 1. ;
106  cIt = electron.superCluster()->clustersBegin() ;
107  cIt != electron.superCluster()->clustersEnd() ;
108  ++cIt )
109  {
110  const reco::CaloClusterPtr cc = *cIt ;
111  crackcor *=
112  ( electron.superCluster()->rawEnergy()
113  + cc->energy()*(crackCorrectionFunction_->getValue(*cc)-1.) )
114  / electron.superCluster()->rawEnergy() ;
115  }
116  newEnergy *= crackcor ;
117 
118  // register final value
119  electron.setCorrectedEcalEnergy(newEnergy) ;
120 
121  }
bool isEcalEnergyCorrected() const
Definition: GsfElectron.h:719
bool isEE() const
Definition: GsfElectron.h:331
bool isEB() const
Definition: GsfElectron.h:330
float fBremEta(float sigmaPhiSigmaEta, float eta, int algorithm, reco::GsfElectron::Classification cl) const
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:168
float fEt(float et, int algorithm, reco::GsfElectron::Classification cl) const
Classification classification() const
Definition: GsfElectron.h:642
return(e1-e2)*(e1-e2)+dp *dp
float fEnergy(float e, int algorithm, reco::GsfElectron::Classification cl) const
void setCorrectedEcalEnergy(float newEnergy)
Definition: GsfElectron.cc:161
float fEta(float energy, float eta, int algorithm) const
EcalClusterFunctionBaseClass * crackCorrectionFunction_
void ElectronEnergyCorrector::classBasedParameterizationUncertainty ( reco::GsfElectron electron)

Definition at line 20 of file ElectronEnergyCorrector.cc.

References reco::GsfElectron::classification(), EnergyUncertaintyElectronSpecific::computeElectronEnergyUncertainty(), reco::GsfElectron::correctedEcalEnergy(), energyError(), reco::GsfElectron::setCorrectedEcalEnergyError(), and reco::GsfElectron::superCluster().

Referenced by GsfElectronAlgo::createElectron().

21  {
23  double energyError = 999. ;
24  double ecalEnergy = electron.correctedEcalEnergy() ;
25  double eleEta = electron.superCluster()->eta() ;
26  double brem = electron.superCluster()->etaWidth()/electron.superCluster()->phiWidth() ;
27  energyError = uncertainty.computeElectronEnergyUncertainty(electron.classification(),eleEta,brem,ecalEnergy) ;
28  electron.setCorrectedEcalEnergyError(energyError) ;
29  }
void setCorrectedEcalEnergyError(float newEnergyError)
Definition: GsfElectron.cc:158
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:168
float correctedEcalEnergy() const
Definition: GsfElectron.h:720
Classification classification() const
Definition: GsfElectron.h:642
double computeElectronEnergyUncertainty(reco::GsfElectron::Classification c, double eta, double brem, double energy)
float energyError(float E, float *par)
float ElectronEnergyCorrector::fBremEta ( float  sigmaPhiSigmaEta,
float  eta,
int  algorithm,
reco::GsfElectron::Classification  cl 
) const
private

Definition at line 159 of file ElectronEnergyCorrector.cc.

References GetRecoTauVFromDQM_MC_cff::cl, reco::GsfElectron::GAP, and fitWZ::par0.

160  {
161  // corrections for electrons
162  if (algorithm!=0)
163  {
164  edm::LogWarning("ElectronEnergyCorrector::fBremEta")<<"algorithm should be 0 for electrons !" ;
165  return 1. ;
166  }
167 
168  const float etaCrackMin = 1.44 ;
169  const float etaCrackMax = 1.56 ;
170 
171  //STD
172  const int nBinsEta = 14 ;
173  float leftEta[nBinsEta] = { 0.02, 0.25, 0.46, 0.81, 0.91, 1.01, 1.16, etaCrackMax, 1.653, 1.8, 2.0, 2.2, 2.3, 2.4 } ;
174  float rightEta[nBinsEta] = { 0.25, 0.42, 0.77, 0.91, 1.01, 1.13, etaCrackMin, 1.653, 1.8, 2.0, 2.2, 2.3, 2.4, 2.5 } ;
175 
176  // eta = 0
177  if ( TMath::Abs(eta) < leftEta[0] ) { eta = 0.02 ; }
178 
179  // outside acceptance
180  if ( TMath::Abs(eta) >= rightEta[nBinsEta-1] ) { eta = 2.49 ; } //if (DBG) std::cout << " WARNING [applyScCorrections]: TMath::Abs(eta) >= rightEta[nBinsEta-1] " << std::endl;}
181 
182  int tmpEta = -1 ;
183  for (int iEta = 0; iEta < nBinsEta; ++iEta)
184  {
185  if ( leftEta[iEta] <= TMath::Abs(eta) && TMath::Abs(eta) <rightEta[iEta] )
186  { tmpEta = iEta ; }
187  }
188 
189  float xcorr[nBinsEta][reco::GsfElectron::GAP+1] =
190  {
191  { 1.00227, 1.00227, 1.00227, 1.00227, 1.00227 },
192  { 1.00252, 1.00252, 1.00252, 1.00252, 1.00252 },
193  { 1.00225, 1.00225, 1.00225, 1.00225, 1.00225 },
194  { 1.00159, 1.00159, 1.00159, 1.00159, 1.00159 },
195  { 0.999475, 0.999475, 0.999475, 0.999475, 0.999475 },
196  { 0.997203, 0.997203, 0.997203, 0.997203, 0.997203 },
197  { 0.993886, 0.993886, 0.993886, 0.993886, 0.993886 },
198  { 0.971262, 0.971262, 0.971262, 0.971262, 0.971262 },
199  { 0.975922, 0.975922, 0.975922, 0.975922, 0.975922 },
200  { 0.979087, 0.979087, 0.979087, 0.979087, 0.979087 },
201  { 0.98495, 0.98495, 0.98495, 0.98495, 0.98495 },
202  { 0.98781, 0.98781, 0.98781, 0.98781, 0.98781 },
203  { 0.989546, 0.989546, 0.989546, 0.989546, 0.989546 },
204  { 0.989638, 0.989638, 0.989638, 0.989638, 0.989638 }
205  } ;
206 
207  float par0[nBinsEta][reco::GsfElectron::GAP+1] =
208  {
209  { 0.994949, 1.00718, 1.00718, 1.00556, 1.00718 },
210  { 1.009, 1.00713, 1.00713, 1.00248, 1.00713 },
211  { 0.999395, 1.00641, 1.00641, 1.00293, 1.00641 },
212  { 0.988662, 1.00761, 1.001, 0.99972, 1.00761 },
213  { 0.998443, 1.00682, 1.001, 1.00282, 1.00682 },
214  { 1.00285, 1.0073, 1.001, 1.00396, 1.0073 },
215  { 0.993053, 1.00462, 1.01341, 1.00184, 1.00462 },
216  { 1.10561, 0.972798, 1.02835, 0.995218, 0.972798 },
217  { 0.893741, 0.981672, 0.98982, 1.01712, 0.981672 },
218  { 0.911123, 0.98251, 1.03466, 1.00824, 0.98251 },
219  { 0.981931, 0.986123, 0.954295, 1.0202, 0.986123 },
220  { 0.905634, 0.990124, 0.928934, 0.998492, 0.990124 },
221  { 0.919343, 0.990187, 0.967526, 0.963923, 0.990187 },
222  { 0.844783, 0.99372, 0.923808, 0.953001, 0.99372 }
223  } ;
224 
225  float par1[nBinsEta][reco::GsfElectron::GAP+1] =
226  {
227  { 0.0111034, -0.00187886, -0.00187886, -0.00289304, -0.00187886 },
228  { -0.00969012, -0.00227574, -0.00227574, -0.00182187, -0.00227574 },
229  { 0.00389454, -0.00259935, -0.00259935, -0.00211059, -0.00259935 },
230  { 0.017095, -0.00433692, -0.00302335, -0.00241385, -0.00433692 },
231  { -0.00049009, -0.00551324, -0.00302335, -0.00532352, -0.00551324 },
232  { -0.00252723, -0.00799669, -0.00302335, -0.00823109, -0.00799669 },
233  { 0.00332567, -0.00870057, -0.0170581, -0.011482, -0.00870057 },
234  { -0.213285, -0.000771577, -0.036007, -0.0187526, -0.000771577 },
235  { 0.1741, -0.00202028, -0.0233995, -0.0302066, -0.00202028 },
236  { 0.152794, 0.00441308, -0.0468563, -0.0158817, 0.00441308 },
237  { 0.0351465, 0.00832913, 0.0358028, -0.0233262, 0.00832913 },
238  { 0.185781, 0.00742879, 0.08858, 0.00568078, 0.00742879 },
239  { 0.153088, 0.0094608, 0.0489979, 0.0491897, 0.0094608 },
240  { 0.296681, 0.00560406, 0.106492, 0.0652007, 0.00560406 }
241  } ;
242 
243  float par2[nBinsEta][reco::GsfElectron::GAP+1] =
244  {
245  { -0.00330844, 0, 0, 5.62441e-05, 0 },
246  { 0.00329373, 0, 0, -0.000113883, 0 },
247  { -0.00104661, 0, 0, -0.000152794, 0 },
248  { -0.0060409, 0, -0.000257724, -0.000202099, 0 },
249  { -0.000742866, 0, -0.000257724, -2.06003e-05, 0 },
250  { -0.00205425, 0, -0.000257724, 3.84179e-05, 0 },
251  { -0.00350757, 0, 0.000590483, 0.000323723, 0 },
252  { 0.0794596, -0.00276696, 0.00205854, 0.000356716, -0.00276696 },
253  { -0.092436, -0.00471028, 0.00062096, 0.00088347, -0.00471028 },
254  { -0.0855029, -0.00809139, 0.00284102, -0.00366903, -0.00809139 },
255  { -0.0306209, -0.00944584, -0.0145892, -0.00176969, -0.00944584 },
256  { -0.0996414, -0.00960462, -0.0328264, -0.00983844, -0.00960462 },
257  { -0.0784107, -0.010172, -0.0256722, -0.0215133, -0.010172 },
258  { -0.145815, -0.00943169, -0.0414525, -0.027087, -0.00943169 }
259  } ;
260 
261  float sigmaPhiSigmaEtaMin[reco::GsfElectron::GAP+1] = { 0.8, 0.8, 0.8, 0.8, 0.8 } ;
262  float sigmaPhiSigmaEtaMax[reco::GsfElectron::GAP+1] = { 5., 5., 5., 5., 5. } ;
263  float sigmaPhiSigmaEtaFit[reco::GsfElectron::GAP+1] = { 1.2, 1.2, 1.2, 1.2, 1.2 } ;
264 
265  // extra protections
266  // fix sigmaPhiSigmaEta boundaries
267  if ( sigmaPhiSigmaEta < sigmaPhiSigmaEtaMin[cl] )
268  { sigmaPhiSigmaEta = sigmaPhiSigmaEtaMin[cl] ; }
269  if ( sigmaPhiSigmaEta > sigmaPhiSigmaEtaMax[cl] )
270  { sigmaPhiSigmaEta = sigmaPhiSigmaEtaMax[cl] ; }
271 
272  // In eta cracks/gaps
273  if ( tmpEta == -1 ) // need to interpolate
274  {
275  float tmpInter = 1 ;
276  for ( int iEta = 0 ; iEta < nBinsEta-1 ; ++iEta )
277  {
278  if ( rightEta[iEta] <= TMath::Abs(eta) && TMath::Abs(eta) <leftEta[iEta+1] )
279  {
280  if ( sigmaPhiSigmaEta >= sigmaPhiSigmaEtaFit[cl] )
281  {
282  tmpInter =
283  ( par0[iEta][cl] +
284  sigmaPhiSigmaEta*par1[iEta][cl] +
285  sigmaPhiSigmaEta*sigmaPhiSigmaEta*par2[iEta][cl] +
286  par0[iEta+1][cl] +
287  sigmaPhiSigmaEta*par1[iEta+1][cl] +
288  sigmaPhiSigmaEta*sigmaPhiSigmaEta*par2[iEta+1][cl] ) / 2. ;
289  }
290  else tmpInter = (xcorr[iEta][cl] + xcorr[iEta+1][cl])/2. ;
291  }
292  }
293  return tmpInter ;
294  }
295 
296  if (sigmaPhiSigmaEta >= sigmaPhiSigmaEtaFit[cl])
297  { return par0[tmpEta][cl] + sigmaPhiSigmaEta*par1[tmpEta][cl] + sigmaPhiSigmaEta*sigmaPhiSigmaEta*par2[tmpEta][cl] ; }
298  else
299  { return xcorr[tmpEta][cl] ; }
300 
301  return 1. ;
302  }
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
T eta() const
tuple par0
Definition: fitWZ.py:51
float ElectronEnergyCorrector::fEnergy ( float  e,
int  algorithm,
reco::GsfElectron::Classification  cl 
) const
private

Definition at line 343 of file ElectronEnergyCorrector.cc.

References GetRecoTauVFromDQM_MC_cff::cl, create_public_lumi_plots::exp, reco::GsfElectron::GAP, and fitWZ::par0.

344  {
345  if (algorithm==0) // Electrons EB
346  { return 1. ; }
347  else if (algorithm==1) // Electrons EE
348  {
349  float par0[reco::GsfElectron::GAP+1] = { 400, 400, 400, 400, 400 } ;
350  float par1[reco::GsfElectron::GAP+1] = { 0.999545, 0.982475, 0.986217, 0.996763, 0.982475 } ;
351  float par2[reco::GsfElectron::GAP+1] = { 1.26568e-05, 4.95413e-05, 5.02161e-05, 2.8485e-06, 4.95413e-05 } ;
352  float par3[reco::GsfElectron::GAP+1] = { 0.0696757, 0.16886, 0.115317, 0.12098, 0.16886 } ;
353  float par4[reco::GsfElectron::GAP+1] = { -54.3468, -30.1517, -26.3805, -62.0538, -30.1517 } ;
354 
355  if ( E > par0[cl] ) { E = par0[cl] ; }
356  if ( E < 0 ) { return 1. ; }
357  if ( 0 <= E && E <= par0[cl] ) { return (par1[cl] + E*par2[cl] )*(1- par3[cl]*exp(E/par4[cl] )) ; }
358  }
359  else
360  { edm::LogWarning("ElectronEnergyCorrector::fEnergy")<<"algorithm should be 0 or 1 for electrons !" ; }
361  return 1.;
362  }
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
tuple par0
Definition: fitWZ.py:51
float ElectronEnergyCorrector::fEt ( float  et,
int  algorithm,
reco::GsfElectron::Classification  cl 
) const
private

Definition at line 304 of file ElectronEnergyCorrector.cc.

References GetRecoTauVFromDQM_MC_cff::cl, create_public_lumi_plots::exp, and reco::GsfElectron::GAP.

305  {
306  if (algorithm==0) //Electrons EB
307  {
308  const float parClassIndep[5] = { 0.97213, 0.999528, 5.61192e-06, 0.0143269, -17.1776 } ;
309  const float par[reco::GsfElectron::GAP+1][5] =
310  {
311  { 0.974327, 0.996127, 5.99401e-05, 0.159813, -3.80392 },
312  { 0.97213, 0.999528, 5.61192e-06, 0.0143269, -17.1776 },
313  { 0.940666, 0.988894, 0.00017474, 0.25603, -4.58153 },
314  { 0.969526, 0.98572, 0.000193842, 4.21548, -1.37159 },
315  { parClassIndep[0], parClassIndep[1], parClassIndep[2], parClassIndep[3], parClassIndep[4] }
316  } ;
317  if ( ET > 200 ) { ET = 200 ; }
318  if ( ET > 100 ) { return (parClassIndep[1]+ET*parClassIndep[2]) * (1-parClassIndep[3]*exp(ET/parClassIndep[4])) ; }
319  if ( ET >= 10 ) { return (par[cl][1]+ET*par[cl][2]) * (1-par[cl][3]*exp(ET/par[cl][4])) ; }
320  if ( ET >=5 ) { return par[cl][0] ; }
321  return 1. ;
322  }
323  else if (algorithm==1) //Electrons EE
324  {
325  float par[reco::GsfElectron::GAP+1][5] =
326  {
327  { 0.930081, 0.996683, 3.54079e-05, 0.0460187, -23.2461 },
328  { 0.930081, 0.996683, 3.54079e-05, 0.0460187, -23.2461 },
329  { 0.930081, 0.996683, 3.54079e-05, 0.0460187, -23.2461 },
330  { 0.930081, 0.996683, 3.54079e-05, 0.0460187, -23.2461 },
331  { 0.930081, 0.996683, 3.54079e-05, 0.0460187, -23.2461 }
332  };
333  if ( ET > 200 ) { ET = 200 ; }
334  if ( ET < 5 ) { return 1. ; }
335  if ( 5 <= ET && ET < 10 ) { return par[cl][0] ; }
336  if ( 10 <= ET && ET <= 200 ) { return ( par[cl][1] + ET*par[cl][2])*(1-par[cl][3]*exp(ET/par[cl][4])) ; }
337  }
338  else
339  { edm::LogWarning("ElectronEnergyCorrector::fEt")<<"algorithm should be 0 or 1 for electrons !" ; }
340  return 1. ;
341  }
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
#define ET
float ElectronEnergyCorrector::fEta ( float  energy,
float  eta,
int  algorithm 
) const
private

Definition at line 128 of file ElectronEnergyCorrector.cc.

References relval_parameters_module::energy, and p1.

129 {
130 
131  // corrections for electrons
132  if (algorithm!=0) {
133  edm::LogWarning("ElectronEnergyCorrector::fEta")<<"algorithm should be 0 for electrons !" ;
134  return energy;
135  }
136  //std::cout << "fEta function" << std::endl;
137 
138  // this correction is setup only for EB
139  if ( algorithm != 0 ) return energy;
140 
141  float ieta = fabs(eta)*(5/0.087);
142  // bypass the DB reading for the time being
143  //float p0 = (params_->params())[0]; // should be 40.2198
144  //float p1 = (params_->params())[1]; // should be -3.03103e-6
145  float p0 = 40.2198 ;
146  float p1 = -3.03103e-6 ;
147 
148  //std::cout << "ieta=" << ieta << std::endl;
149 
150  float correctedEnergy = energy;
151  if ( ieta < p0 ) correctedEnergy = energy;
152  else correctedEnergy = energy/(1.0 + p1*(ieta-p0)*(ieta-p0));
153  //std::cout << "ECEC fEta = " << correctedEnergy << std::endl;
154  return correctedEnergy;
155 
156  }
< trclass="colgroup">< tdclass="colgroup"colspan=5 > Ecal cluster collections</td ></tr >< tr >< td >< ahref="classreco_1_1BasicCluster.html"> reco::BasicCluster</a ></td >< td >< ahref="DataFormats_EgammaReco.html"> reco::BasicClusterCollection</a ></td >< td >< ahref="#"> hybridSuperClusters</a ></td >< tdclass="description"> Basic clusters reconstructed with hybrid algorithm(barrel only)</td >< td >S.Rahatlou</td ></tr >< tr >< td >< a href
T eta() const
double p1[4]
Definition: TauolaWrapper.h:89
double ElectronEnergyCorrector::fEtaBarrelBad ( double  scEta) const
private

Definition at line 382 of file ElectronEnergyCorrector.cc.

References funct::abs(), p1, p2, p3, p4, and x.

383  {
384  // f(eta) for the class = 30 (estimated from 1Mevt single e sample)
385  // Ivica's new corrections 01/06
386  float p0 = 9.99063e-01;
387  float p1 = -2.63341e-02;
388  float p2 = 5.16054e-02;
389  float p3 = -4.95976e-02;
390  float p4 = 3.62304e-03;
391  double x = (double) std::abs(scEta) ;
392  return p0 + p1*x + p2*x*x + p3*x*x*x + p4*x*x*x*x ;
393  }
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
double p1[4]
Definition: TauolaWrapper.h:89
Definition: DDAxes.h:10
double p3[4]
Definition: TauolaWrapper.h:91
double ElectronEnergyCorrector::fEtaBarrelGood ( double  scEta) const
private

Definition at line 369 of file ElectronEnergyCorrector.cc.

References funct::abs(), p1, p2, p3, p4, and x.

370  {
371  // f(eta) for the first 3 classes (0, 10 and 20) (estimated from 1Mevt single e sample)
372  // Ivica's new corrections 01/06
373  float p0 = 1.00149e+00 ;
374  float p1 = -2.06622e-03 ;
375  float p2 = -1.08793e-02 ;
376  float p3 = 1.54392e-02 ;
377  float p4 = -1.02056e-02 ;
378  double x = (double) std::abs(scEta) ;
379  return p0 + p1*x + p2*x*x + p3*x*x*x + p4*x*x*x*x ;
380  }
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
double p1[4]
Definition: TauolaWrapper.h:89
Definition: DDAxes.h:10
double p3[4]
Definition: TauolaWrapper.h:91
double ElectronEnergyCorrector::fEtaEndcapBad ( double  scEta) const
private

Definition at line 408 of file ElectronEnergyCorrector.cc.

References funct::abs(), p1, p2, p3, p4, and x.

409  {
410  // f(eta) for the class = 130-134
411  // Ivica's new corrections 01/06
412  float p0 = -4.25221e+00 ;
413  float p1 = 1.01936e+01 ;
414  float p2 = -7.48247e+00 ;
415  float p3 = 2.45520e+00 ;
416  float p4 = -3.02872e-01 ;
417  double x = (double) std::abs(scEta);
418  return p0 + p1*x + p2*x*x + p3*x*x*x + p4*x*x*x*x ;
419  }
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
double p1[4]
Definition: TauolaWrapper.h:89
Definition: DDAxes.h:10
double p3[4]
Definition: TauolaWrapper.h:91
double ElectronEnergyCorrector::fEtaEndcapGood ( double  scEta) const
private

Definition at line 395 of file ElectronEnergyCorrector.cc.

References funct::abs(), p1, p2, p3, p4, and x.

396  {
397  // f(eta) for the first 3 classes (100, 110 and 120)
398  // Ivica's new corrections 01/06
399  float p0 = -8.51093e-01 ;
400  float p1 = 3.54266e+00 ;
401  float p2 = -2.59288e+00 ;
402  float p3 = 8.58945e-01 ;
403  float p4 = -1.07844e-01 ;
404  double x = (double) std::abs(scEta) ;
405  return p0 + p1*x + p2*x*x + p3*x*x*x + p4*x*x*x*x ;
406  }
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
double p1[4]
Definition: TauolaWrapper.h:89
Definition: DDAxes.h:10
double p3[4]
Definition: TauolaWrapper.h:91
void ElectronEnergyCorrector::simpleParameterizationUncertainty ( reco::GsfElectron electron)

Definition at line 34 of file ElectronEnergyCorrector.cc.

References reco::GsfElectron::correctedEcalEnergy(), energyError(), relativeConstraints::error, reco::GsfElectron::isEB(), reco::GsfElectron::isEE(), and reco::GsfElectron::setCorrectedEcalEnergyError().

Referenced by GsfElectronAlgo::createElectron().

35  {
36  double error = 999. ;
37  double ecalEnergy = electron.correctedEcalEnergy() ;
38 
39  if (electron.isEB())
40  {
41  float parEB[3] = { 5.24e-02, 2.01e-01, 1.00e-02} ;
42  error = ecalEnergy*energyError(ecalEnergy,parEB) ;
43  }
44  else if (electron.isEE())
45  {
46  float parEE[3] = { 1.46e-01, 9.21e-01, 1.94e-03} ;
47  error = ecalEnergy*energyError(ecalEnergy,parEE) ;
48  }
49  else
50  { edm::LogWarning("ElectronEnergyCorrector::simpleParameterizationUncertainty")<<"nor barrel neither endcap electron !" ; }
51 
52  electron.setCorrectedEcalEnergyError(error) ;
53  }
bool isEE() const
Definition: GsfElectron.h:331
bool isEB() const
Definition: GsfElectron.h:330
void setCorrectedEcalEnergyError(float newEnergyError)
Definition: GsfElectron.cc:158
float correctedEcalEnergy() const
Definition: GsfElectron.h:720
float energyError(float E, float *par)

Member Data Documentation

EcalClusterFunctionBaseClass* ElectronEnergyCorrector::crackCorrectionFunction_
private

Definition at line 36 of file ElectronEnergyCorrector.h.