CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EcalClusterEnergyCorrection Class Reference
Inheritance diagram for EcalClusterEnergyCorrection:
EcalClusterFunctionBaseClass

Public Member Functions

void checkInit () const
 
 EcalClusterEnergyCorrection (const edm::ParameterSet &, edm::ConsumesCollector iC)
 
const EcalClusterEnergyCorrectionParametersgetParameters () const
 
float getValue (const reco::SuperCluster &, const int mode) const override
 
float getValue (const reco::BasicCluster &, const EcalRecHitCollection &) const override
 
void init (const edm::EventSetup &es) override
 
- Public Member Functions inherited from EcalClusterFunctionBaseClass
virtual float getValue (const reco::CaloCluster &) const
 
virtual ~EcalClusterFunctionBaseClass ()
 

Private Member Functions

float fBrem (float e, float eta, int algorithm) const
 
float fEta (float e, float eta, int algorithm) const
 
float fEtEta (float et, float eta, int algorithm) const
 

Private Attributes

const EcalClusterEnergyCorrectionParametersparams_ = nullptr
 
const edm::ESGetToken< EcalClusterEnergyCorrectionParameters, EcalClusterEnergyCorrectionParametersRcdparamsToken_
 

Detailed Description

Function that provides supercluster energy correction due to Bremsstrahlung loss

$Id: EcalClusterEnergyCorrection.h $Date: $Revision:

Author
Yurii Maravin, KSU, March 2009

Definition at line 17 of file EcalClusterEnergyCorrection.cc.

Constructor & Destructor Documentation

◆ EcalClusterEnergyCorrection()

EcalClusterEnergyCorrection::EcalClusterEnergyCorrection ( const edm::ParameterSet ,
edm::ConsumesCollector  iC 
)
inline

Definition at line 19 of file EcalClusterEnergyCorrection.cc.

19 : paramsToken_(iC.esConsumes()) {}
const edm::ESGetToken< EcalClusterEnergyCorrectionParameters, EcalClusterEnergyCorrectionParametersRcd > paramsToken_

Member Function Documentation

◆ checkInit()

void EcalClusterEnergyCorrection::checkInit ( ) const

Definition at line 246 of file EcalClusterEnergyCorrection.cc.

References Exception, and params_.

Referenced by getValue().

246  {
247  if (!params_) {
248  // non-initialized function parameters: throw exception
249  throw cms::Exception("EcalClusterEnergyCorrection::checkInit()")
250  << "Trying to access an uninitialized crack correction function.\n"
251  "Please call `init( edm::EventSetup &)' before any use of the function.\n";
252  }
253 }
const EcalClusterEnergyCorrectionParameters * params_

◆ fBrem()

float EcalClusterEnergyCorrection::fBrem ( float  e,
float  eta,
int  algorithm 
) const
private

Definition at line 62 of file EcalClusterEnergyCorrection.cc.

References a, qcdUeDQM_cfi::algorithm, b, DummyCfis::c, MillePedeFileConverter_cfg::e, hltrates_dqm_sourceclient-live_cfg::offset, LaserDQM_cfg::p1, SiStripOfflineCRack_cfg::p2, chargedHadronTrackResolutionFilter_cfi::p3, EcalFunParams::params(), params_, DiMuonV_cfg::threshold, and y.

Referenced by getValue().

62  {
63  // brem == phiWidth/etaWidth of the SuperCluster
64  // e == energy of the SuperCluster
65  // first parabola (for br > threshold)
66  // p0 + p1*x + p2*x^2
67  // second parabola (for br <= threshold)
68  // ax^2 + bx + c, make y and y' the same in threshold
69  // y = p0 + p1*threshold + p2*threshold^2
70  // yprime = p1 + 2*p2*threshold
71  // a = p3
72  // b = yprime - 2*a*threshold
73  // c = y - a*threshold^2 - b*threshold
74 
75  int offset;
76  if (algorithm == 0)
77  offset = 0;
78  else if (algorithm == 1)
79  offset = 20;
80  else {
81  // not supported, produce no correction
82  return e;
83  }
84 
85  //Make No Corrections if brem is invalid!
86  if (brem == 0)
87  return e;
88 
89  float bremLowThr = (params_->params())[2 + offset];
90  float bremHighThr = (params_->params())[3 + offset];
91  if (brem < bremLowThr)
92  brem = bremLowThr;
93  if (brem > bremHighThr)
94  brem = bremHighThr;
95 
96  // Parameters provided in cfg file
97  float p0 = (params_->params())[4 + offset];
98  float p1 = (params_->params())[5 + offset];
99  float p2 = (params_->params())[6 + offset];
100  float p3 = (params_->params())[7 + offset];
101  float p4 = (params_->params())[8 + offset];
102  //
103  float threshold = p4;
104 
105  float y = p0 * threshold * threshold + p1 * threshold + p2;
106  float yprime = 2 * p0 * threshold + p1;
107  float a = p3;
108  float b = yprime - 2 * a * threshold;
109  float c = y - a * threshold * threshold - b * threshold;
110 
111  float fCorr = 1;
112  if (brem < threshold)
113  fCorr = p0 * brem * brem + p1 * brem + p2;
114  else
115  fCorr = a * brem * brem + b * brem + c;
116 
117  //std::cout << "ECEC fBrem " << e/fCorr << std::endl;
118  return e / fCorr;
119 }
const EcalClusterEnergyCorrectionParameters * params_
EcalFunctionParameters & params()
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121

◆ fEta()

float EcalClusterEnergyCorrection::fEta ( float  e,
float  eta,
int  algorithm 
) const
private

Definition at line 44 of file EcalClusterEnergyCorrection.cc.

References qcdUeDQM_cfi::algorithm, hcalRecHitTable_cff::energy, PVValHelper::eta, hcalRecHitTable_cff::ieta, LaserDQM_cfg::p1, EcalFunParams::params(), and params_.

Referenced by getValue().

44  {
45  // this correction is setup only for EB
46  if (algorithm != 0)
47  return energy;
48 
49  float ieta = fabs(eta) * (5 / 0.087);
50  float p0 = (params_->params())[0]; // should be 40.2198
51  float p1 = (params_->params())[1]; // should be -3.03103e-6
52 
53  float correctedEnergy = energy;
54  if (ieta < p0)
55  correctedEnergy = energy;
56  else
57  correctedEnergy = energy / (1.0 + p1 * (ieta - p0) * (ieta - p0));
58  //std::cout << "ECEC fEta = " << correctedEnergy << std::endl;
59  return correctedEnergy;
60 }
const EcalClusterEnergyCorrectionParameters * params_
EcalFunctionParameters & params()

◆ fEtEta()

float EcalClusterEnergyCorrection::fEtEta ( float  et,
float  eta,
int  algorithm 
) const
private

Definition at line 121 of file EcalClusterEnergyCorrection.cc.

References qcdUeDQM_cfi::algorithm, l1tnanotables_cff::et, PVValHelper::eta, hltrates_dqm_sourceclient-live_cfg::offset, LaserDQM_cfg::p1, SiStripOfflineCRack_cfg::p2, chargedHadronTrackResolutionFilter_cfi::p3, EcalFunParams::params(), params_, and mathSSE::sqrt().

Referenced by getValue().

121  {
122  // et -- Et of the SuperCluster (with respect to (0,0,0))
123  // eta -- eta of the SuperCluster
124 
125  //std::cout << "fEtEta, mode = " << algorithm << std::endl;
126  //std::cout << "ECEC: p0 " << (params_->params())[9] << " " << (params_->params())[10] << " " << (params_->params())[11] << " " << (params_->params())[12] << std::endl;
127  //std::cout << "ECEC: p1 " << (params_->params())[13] << " " << (params_->params())[14] << " " << (params_->params())[15] << " " << (params_->params())[16] << std::endl;
128  //std::cout << "ECEC: fcorr " << (params_->params())[17] << " " << (params_->params())[18] << " " << (params_->params())[19] << std::endl;
129 
130  float fCorr = 0.;
131  int offset;
132  if (algorithm == 0)
133  offset = 0;
134  else if (algorithm == 1)
135  offset = 20;
136  else if (algorithm == 10)
137  offset = 28;
138  else if (algorithm == 11)
139  offset = 39;
140  else {
141  // not supported, produce no correction
142  return et;
143  }
144 
145  // Barrel
146  if (algorithm == 0 || algorithm == 10) {
147  float p0 = (params_->params())[9 + offset] +
148  (params_->params())[10 + offset] / (et + (params_->params())[11 + offset]) +
149  (params_->params())[12 + offset] / (et * et);
150  float p1 = (params_->params())[13 + offset] +
151  (params_->params())[14 + offset] / (et + (params_->params())[15 + offset]) +
152  (params_->params())[16 + offset] / (et * et);
153 
154  fCorr = p0 + p1 * atan((params_->params())[17 + offset] * ((params_->params())[18 + offset] - fabs(eta))) +
155  (params_->params())[19 + offset] * fabs(eta);
156 
157  } else if (algorithm == 1 || algorithm == 11) { // Endcap
158  float p0 = (params_->params())[9 + offset] + (params_->params())[10 + offset] / sqrt(et);
159  float p1 = (params_->params())[11 + offset] + (params_->params())[12 + offset] / sqrt(et);
160  float p2 = (params_->params())[13 + offset] + (params_->params())[14 + offset] / sqrt(et);
161  float p3 = (params_->params())[15 + offset] + (params_->params())[16 + offset] / sqrt(et);
162 
163  fCorr = p0 + p1 * fabs(eta) + p2 * eta * eta + p3 / fabs(eta);
164  }
165 
166  // cap the correction at 50%
167  if (fCorr < 0.5)
168  fCorr = 0.5;
169  if (fCorr > 1.5)
170  fCorr = 1.5;
171 
172  //std::cout << "ECEC fEtEta " << et/fCorr << std::endl;
173  return et / fCorr;
174 }
const EcalClusterEnergyCorrectionParameters * params_
T sqrt(T t)
Definition: SSEVec.h:23
EcalFunctionParameters & params()

◆ getParameters()

const EcalClusterEnergyCorrectionParameters* EcalClusterEnergyCorrection::getParameters ( ) const
inline

Definition at line 22 of file EcalClusterEnergyCorrection.cc.

References params_.

22 { return params_; }
const EcalClusterEnergyCorrectionParameters * params_

◆ getValue() [1/2]

float EcalClusterEnergyCorrection::getValue ( const reco::SuperCluster superCluster,
const int  mode 
) const
overridevirtual

Implements EcalClusterFunctionBaseClass.

Definition at line 176 of file EcalClusterEnergyCorrection.cc.

References qcdUeDQM_cfi::algorithm, checkInit(), hcalRecHitTable_cff::energy, reco::CaloCluster::energy(), PVValHelper::eta, reco::CaloCluster::eta(), reco::SuperCluster::etaWidth(), fBrem(), fEta(), fEtEta(), ALCARECOPromptCalibProdSiPixelAli0T_cff::mode, reco::SuperCluster::phiWidth(), reco::SuperCluster::preshowerEnergy(), and reco::SuperCluster::rawEnergy().

176  {
177  // mode flags:
178  // hybrid modes: 1 - return f(eta) correction in GeV
179  // 2 - return f(eta) + f(brem) correction
180  // 3 - return f(eta) + f(brem) + f(et, eta) correction
181  // multi5x5: 4 - return f(brem) correction
182  // 5 - return f(brem) + f(et, eta) correction
183 
184  // special cases: mode = 10 -- return f(et, eta) correction with respect to already corrected SC in barrel
185  // mode = 11 -- return f(et, eta) correction with respect to already corrected SC in endcap
186 
187  checkInit();
188 
189  float eta = fabs(superCluster.eta());
190  float brem = superCluster.phiWidth() / superCluster.etaWidth();
191  int algorithm = -1;
192 
193  if (mode <= 3 || mode == 10) {
194  // algorithm: hybrid
195  algorithm = 0;
196 
197  float energy = superCluster.rawEnergy();
198  if (mode == 10) {
199  algorithm = 10;
200  energy = superCluster.energy();
201  }
202  float correctedEnergy = fEta(energy, eta, algorithm);
203 
204  if (mode == 1) {
205  return correctedEnergy - energy;
206 
207  } else {
208  // now apply F(brem)
209  correctedEnergy = fBrem(correctedEnergy, brem, algorithm);
210  if (mode == 2) {
211  return correctedEnergy - energy;
212  }
213 
214  float correctedEt = correctedEnergy / cosh(eta);
215  correctedEt = fEtEta(correctedEt, eta, algorithm);
216  correctedEnergy = correctedEt * cosh(eta);
217  return correctedEnergy - energy;
218  }
219  } else if (mode == 4 || mode == 5 || mode == 11) {
220  algorithm = 1;
221 
222  float energy = superCluster.rawEnergy() + superCluster.preshowerEnergy();
223  if (mode == 11) {
224  algorithm = 11;
225  energy = superCluster.energy();
226  }
227 
228  float correctedEnergy = fBrem(energy, brem, algorithm);
229  if (mode == 4) {
230  return correctedEnergy - energy;
231  }
232 
233  float correctedEt = correctedEnergy / cosh(eta);
234  correctedEt = fEtEta(correctedEt, eta, algorithm);
235  correctedEnergy = correctedEt * cosh(eta);
236  return correctedEnergy - energy;
237 
238  } else {
239  // perform no correction
240  return 0;
241  }
242 }
float fEta(float e, float eta, int algorithm) const
float fEtEta(float et, float eta, int algorithm) const
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:60
double phiWidth() const
obtain phi and eta width of the Super Cluster
Definition: SuperCluster.h:68
float fBrem(float e, float eta, int algorithm) const
double energy() const
cluster energy
Definition: CaloCluster.h:149
double etaWidth() const
Definition: SuperCluster.h:69
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:63

◆ getValue() [2/2]

float EcalClusterEnergyCorrection::getValue ( const reco::BasicCluster ,
const EcalRecHitCollection  
) const
inlineoverridevirtual

Implements EcalClusterFunctionBaseClass.

Definition at line 28 of file EcalClusterEnergyCorrection.cc.

28 { return 0.; };

◆ init()

void EcalClusterEnergyCorrection::init ( const edm::EventSetup es)
overridevirtual

Implements EcalClusterFunctionBaseClass.

Definition at line 244 of file EcalClusterEnergyCorrection.cc.

References edm::EventSetup::getData(), params_, and paramsToken_.

244 { params_ = &es.getData(paramsToken_); }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const EcalClusterEnergyCorrectionParameters * params_
const edm::ESGetToken< EcalClusterEnergyCorrectionParameters, EcalClusterEnergyCorrectionParametersRcd > paramsToken_

Member Data Documentation

◆ params_

const EcalClusterEnergyCorrectionParameters* EcalClusterEnergyCorrection::params_ = nullptr
private

Definition at line 39 of file EcalClusterEnergyCorrection.cc.

Referenced by checkInit(), fBrem(), fEta(), fEtEta(), getParameters(), and init().

◆ paramsToken_

const edm::ESGetToken<EcalClusterEnergyCorrectionParameters, EcalClusterEnergyCorrectionParametersRcd> EcalClusterEnergyCorrection::paramsToken_
private

Definition at line 38 of file EcalClusterEnergyCorrection.cc.

Referenced by init().