CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
PFClusterEMEnergyCorrector Class Reference

#include <PFClusterEMEnergyCorrector.h>

Public Member Functions

void correctEnergies (const edm::Event &evt, const edm::EventSetup &es, const reco::PFCluster::EEtoPSAssociation &assoc, reco::PFClusterCollection &cs)
 
PFClusterEMEnergyCorrectoroperator= (const PFClusterEMEnergyCorrector &)=delete
 
 PFClusterEMEnergyCorrector (const edm::ParameterSet &conf, edm::ConsumesCollector &&cc)
 
 PFClusterEMEnergyCorrector (const PFClusterEMEnergyCorrector &)=delete
 

Private Member Functions

void getAssociatedPSEnergy (const size_t clusIdx, const reco::PFCluster::EEtoPSAssociation &assoc, float &e1, float &e2)
 

Private Attributes

bool applyCrackCorrections_
 
bool applyMVACorrections_
 
bool autoDetectBunchSpacing_
 
edm::EDGetTokenT< unsigned int > bunchSpacing_
 
int bunchSpacingManual_
 
std::unique_ptr< PFEnergyCalibrationcalibrator_
 
std::vector< std::string > condnames_mean_
 
std::vector< std::string > condnames_mean_25ns_
 
std::vector< std::string > condnames_mean_50ns_
 
std::vector< std::string > condnames_sigma_
 
std::vector< std::string > condnames_sigma_25ns_
 
std::vector< std::string > condnames_sigma_50ns_
 
edm::EDGetTokenT< EBSrFlagCollectionebSrFlagToken_
 
edm::EDGetTokenT< EESrFlagCollectioneeSrFlagToken_
 
double maxPtForMVAEvaluation_
 
double meanlimhighEB_
 
double meanlimhighEE_
 
double meanlimlowEB_
 
double meanlimlowEE_
 
double meanoffsetEB_
 
double meanoffsetEE_
 
double meanscaleEB_
 
double meanscaleEE_
 
edm::EDGetTokenT< EcalRecHitCollectionrecHitsEB_
 
edm::EDGetTokenT< EcalRecHitCollectionrecHitsEE_
 
bool setEnergyUncertainty_
 
double sigmalimhighEB_
 
double sigmalimhighEE_
 
double sigmalimlowEB_
 
double sigmalimlowEE_
 
double sigmaoffsetEB_
 
double sigmaoffsetEE_
 
double sigmascaleEB_
 
double sigmascaleEE_
 
bool srfAwareCorrection_
 

Detailed Description

Definition at line 26 of file PFClusterEMEnergyCorrector.h.

Constructor & Destructor Documentation

PFClusterEMEnergyCorrector::PFClusterEMEnergyCorrector ( const edm::ParameterSet conf,
edm::ConsumesCollector &&  cc 
)

Definition at line 14 of file PFClusterEMEnergyCorrector.cc.

References applyCrackCorrections_, applyMVACorrections_, autoDetectBunchSpacing_, bunchSpacing_, bunchSpacingManual_, condnames_mean_, condnames_mean_25ns_, condnames_mean_50ns_, condnames_sigma_, condnames_sigma_25ns_, condnames_sigma_50ns_, ebSrFlagToken_, eeSrFlagToken_, edm::ParameterSet::getParameter(), HLT_2018_cff::InputTag, maxPtForMVAEvaluation_, meanlimhighEB_, meanlimhighEE_, meanlimlowEB_, meanlimlowEE_, meanoffsetEB_, meanoffsetEE_, meanscaleEB_, meanscaleEE_, recHitsEB_, recHitsEE_, setEnergyUncertainty_, sigmalimhighEB_, sigmalimhighEE_, sigmalimlowEB_, sigmalimlowEE_, sigmaoffsetEB_, sigmaoffsetEE_, sigmascaleEB_, sigmascaleEE_, and srfAwareCorrection_.

16  applyCrackCorrections_ = conf.getParameter<bool>("applyCrackCorrections");
17  applyMVACorrections_ = conf.getParameter<bool>("applyMVACorrections");
18  srfAwareCorrection_ = conf.getParameter<bool>("srfAwareCorrection");
19  setEnergyUncertainty_ = conf.getParameter<bool>("setEnergyUncertainty");
20  maxPtForMVAEvaluation_ = conf.getParameter<double>("maxPtForMVAEvaluation");
21 
22  if (applyMVACorrections_) {
23  meanlimlowEB_ = -0.336;
24  meanlimhighEB_ = 0.916;
27 
28  meanlimlowEE_ = -0.336;
29  meanlimhighEE_ = 0.916;
32 
33  sigmalimlowEB_ = 0.001;
34  sigmalimhighEB_ = 0.4;
37 
38  sigmalimlowEE_ = 0.001;
39  sigmalimhighEE_ = 0.4;
42 
43  recHitsEB_ = cc.consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitsEBLabel"));
44  recHitsEE_ = cc.consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitsEELabel"));
45  autoDetectBunchSpacing_ = conf.getParameter<bool>("autoDetectBunchSpacing");
46 
48  bunchSpacing_ = cc.consumes<unsigned int>(edm::InputTag("bunchSpacingProducer"));
50  } else {
51  bunchSpacingManual_ = conf.getParameter<int>("bunchSpacing");
52  }
53 
55  {"ecalPFClusterCorV2_EB_pfSize1_mean_25ns",
56  "ecalPFClusterCorV2_EB_pfSize2_mean_25ns",
57  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_mean_25ns",
58  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_mean_25ns",
59  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_mean_25ns",
60  "ecalPFClusterCorV2_EE_pfSize1_mean_25ns",
61  "ecalPFClusterCorV2_EE_pfSize2_mean_25ns",
62  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_mean_25ns",
63  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_mean_25ns",
64  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_mean_25ns"});
66  {"ecalPFClusterCorV2_EB_pfSize1_sigma_25ns",
67  "ecalPFClusterCorV2_EB_pfSize2_sigma_25ns",
68  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_sigma_25ns",
69  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_sigma_25ns",
70  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_sigma_25ns",
71  "ecalPFClusterCorV2_EE_pfSize1_sigma_25ns",
72  "ecalPFClusterCorV2_EE_pfSize2_sigma_25ns",
73  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_sigma_25ns",
74  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_sigma_25ns",
75  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_sigma_25ns"});
77  {"ecalPFClusterCorV2_EB_pfSize1_mean_50ns",
78  "ecalPFClusterCorV2_EB_pfSize2_mean_50ns",
79  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_mean_50ns",
80  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_mean_50ns",
81  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_mean_50ns",
82  "ecalPFClusterCorV2_EE_pfSize1_mean_50ns",
83  "ecalPFClusterCorV2_EE_pfSize2_mean_50ns",
84  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_mean_50ns",
85  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_mean_50ns",
86  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_mean_50ns"});
88  {"ecalPFClusterCorV2_EB_pfSize1_sigma_50ns",
89  "ecalPFClusterCorV2_EB_pfSize2_sigma_50ns",
90  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_sigma_50ns",
91  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_sigma_50ns",
92  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_sigma_50ns",
93  "ecalPFClusterCorV2_EE_pfSize1_sigma_50ns",
94  "ecalPFClusterCorV2_EE_pfSize2_sigma_50ns",
95  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_sigma_50ns",
96  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_sigma_50ns",
97  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_sigma_50ns"});
98 
99  if (srfAwareCorrection_) {
100  sigmalimlowEE_ = 0.001;
101  sigmalimhighEE_ = 0.1;
104 
107 
108  condnames_mean_.insert(condnames_mean_.end(),
109  {"ecalPFClusterCor2017V2_EB_ZS_mean_25ns",
110  "ecalPFClusterCor2017V2_EB_Full_ptbin1_mean_25ns",
111  "ecalPFClusterCor2017V2_EB_Full_ptbin2_mean_25ns",
112  "ecalPFClusterCor2017V2_EB_Full_ptbin3_mean_25ns",
113  "ecalPFClusterCor2017V2_EE_ZS_mean_25ns",
114  "ecalPFClusterCor2017V2_EE_Full_ptbin1_mean_25ns",
115  "ecalPFClusterCor2017V2_EE_Full_ptbin2_mean_25ns",
116  "ecalPFClusterCor2017V2_EE_Full_ptbin3_mean_25ns"});
117 
118  condnames_sigma_.insert(condnames_sigma_.end(),
119  {"ecalPFClusterCor2017V2_EB_ZS_sigma_25ns",
120  "ecalPFClusterCor2017V2_EB_Full_ptbin1_sigma_25ns",
121  "ecalPFClusterCor2017V2_EB_Full_ptbin2_sigma_25ns",
122  "ecalPFClusterCor2017V2_EB_Full_ptbin3_sigma_25ns",
123  "ecalPFClusterCor2017V2_EE_ZS_sigma_25ns",
124  "ecalPFClusterCor2017V2_EE_Full_ptbin1_sigma_25ns",
125  "ecalPFClusterCor2017V2_EE_Full_ptbin2_sigma_25ns",
126  "ecalPFClusterCor2017V2_EE_Full_ptbin3_sigma_25ns"});
127  }
128  }
129 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< EESrFlagCollection > eeSrFlagToken_
std::unique_ptr< PFEnergyCalibration > calibrator_
edm::EDGetTokenT< unsigned int > bunchSpacing_
edm::EDGetTokenT< EcalRecHitCollection > recHitsEE_
std::vector< std::string > condnames_mean_50ns_
std::vector< std::string > condnames_sigma_25ns_
std::vector< std::string > condnames_mean_25ns_
std::vector< std::string > condnames_sigma_
edm::EDGetTokenT< EcalRecHitCollection > recHitsEB_
edm::EDGetTokenT< EBSrFlagCollection > ebSrFlagToken_
std::vector< std::string > condnames_sigma_50ns_
std::vector< std::string > condnames_mean_
PFClusterEMEnergyCorrector::PFClusterEMEnergyCorrector ( const PFClusterEMEnergyCorrector )
delete

Member Function Documentation

void PFClusterEMEnergyCorrector::correctEnergies ( const edm::Event evt,
const edm::EventSetup es,
const reco::PFCluster::EEtoPSAssociation assoc,
reco::PFClusterCollection cs 
)

a hit can be ZS or forced ZS. A hit can be in Full readout or Forced to be FULL readout if it is ZS, then clusFlag (in binary) = 0001 if it is forced ZS, then clusFlag (in binary) = 0101 if it is FR, then clusFlag (in binary) = 0011 if it is forced FR, then clusFlag (in binary) = 0111 i.e 3rd bit is set. Even if it is forced, we should mark it is as ZS or FR. To take care of it, just check the LSB and second LSB(SLSB)

it is clusFlag==1, 5

Definition at line 131 of file PFClusterEMEnergyCorrector.cc.

References funct::abs(), applyCrackCorrections_, applyMVACorrections_, autoDetectBunchSpacing_, bunchSpacing_, bunchSpacingManual_, calibrator_, condnames_mean_, condnames_mean_25ns_, condnames_mean_50ns_, condnames_sigma_, condnames_sigma_25ns_, condnames_sigma_50ns_, MillePedeFileConverter_cfg::e, ebSrFlagToken_, PFLayer::ECAL_BARREL, eeSrFlagToken_, edm::SortedCollection< T, SORT >::end(), reco::PFCluster::energy(), JetChargeProducer_cfi::exp, edm::SortedCollection< T, SORT >::find(), edm::EventSetup::get(), getAssociatedPSEnergy(), edm::Event::getByToken(), GBRForestD::GetResponse(), training_settings::idx, EBDetId::ieta(), edm::HandleBase::isValid(), EEDetId::ix(), reco::PFCluster::layer(), LogDebug, maxPtForMVAEvaluation_, SiStripPI::mean, meanoffsetEB_, meanoffsetEE_, meanscaleEB_, meanscaleEE_, min(), DiDispStaMuonMonitor_cfi::pt, reco::PFCluster::pt(), EcalReadoutTools::readOutUnitOf(), recHitsEB_, recHitsEE_, reco::CaloCluster::seed(), reco::CaloCluster::setCorrectedEnergy(), reco::CaloCluster::setCorrectedEnergyUncertainty(), setEnergyUncertainty_, sigmaoffsetEB_, sigmaoffsetEE_, sigmascaleEB_, sigmascaleEE_, findQualityFiles::size, and srfAwareCorrection_.

134  {
135  // First deal with pre-MVA corrections
136  // Kept for backward compatibility (and for HLT)
137  if (!applyMVACorrections_) {
138  for (unsigned int idx = 0; idx < cs.size(); ++idx) {
139  reco::PFCluster &cluster = cs[idx];
140  bool iseb = cluster.layer() == PFLayer::ECAL_BARREL;
141  float ePS1 = 0., ePS2 = 0.;
142  if (!iseb)
143  getAssociatedPSEnergy(idx, assoc, ePS1, ePS2);
144  double correctedEnergy = calibrator_->energyEm(cluster, ePS1, ePS2, applyCrackCorrections_);
145  cluster.setCorrectedEnergy(correctedEnergy);
146  }
147  return;
148  }
149 
150  // Common objects for SRF-aware and old style corrections
151  EcalClusterLazyTools lazyTool(evt, es, recHitsEB_, recHitsEE_);
152  EcalReadoutTools readoutTool(evt, es);
153 
154  if (!srfAwareCorrection_) {
155  int bunchspacing = 450;
157  edm::Handle<unsigned int> bunchSpacingH;
158  evt.getByToken(bunchSpacing_, bunchSpacingH);
159  bunchspacing = *bunchSpacingH;
160  } else {
161  bunchspacing = bunchSpacingManual_;
162  }
163 
164  const std::vector<std::string> &condnames_mean = (bunchspacing == 25) ? condnames_mean_25ns_ : condnames_mean_50ns_;
165  const std::vector<std::string> &condnames_sigma =
166  (bunchspacing == 25) ? condnames_sigma_25ns_ : condnames_sigma_50ns_;
167  const unsigned int ncor = condnames_mean.size();
168  std::vector<edm::ESHandle<GBRForestD> > forestH_mean(ncor);
169  std::vector<edm::ESHandle<GBRForestD> > forestH_sigma(ncor);
170 
171  for (unsigned int icor = 0; icor < ncor; ++icor) {
172  es.get<GBRDWrapperRcd>().get(condnames_mean[icor], forestH_mean[icor]);
173  es.get<GBRDWrapperRcd>().get(condnames_sigma[icor], forestH_sigma[icor]);
174  }
175 
176  std::array<float, 5> eval;
177  for (unsigned int idx = 0; idx < cs.size(); ++idx) {
178  reco::PFCluster &cluster = cs[idx];
179  bool iseb = cluster.layer() == PFLayer::ECAL_BARREL;
180  float ePS1 = 0., ePS2 = 0.;
181  if (!iseb)
182  getAssociatedPSEnergy(idx, assoc, ePS1, ePS2);
183 
184  double e = cluster.energy();
185  double pt = cluster.pt();
186  double evale = e;
188  evale *= maxPtForMVAEvaluation_ / pt;
189  }
190  double invE = (e == 0.) ? 0. : 1. / e; //guard against dividing by 0.
191  int size = lazyTool.n5x5(cluster);
192 
193  int ietaix = 0;
194  int iphiiy = 0;
195  if (iseb) {
196  EBDetId ebseed(cluster.seed());
197  ietaix = ebseed.ieta();
198  iphiiy = ebseed.iphi();
199  } else {
200  EEDetId eeseed(cluster.seed());
201  ietaix = eeseed.ix();
202  iphiiy = eeseed.iy();
203  }
204 
205  //find index of corrections (0-4 for EB, 5-9 for EE, depending on cluster size and raw pt)
206  int coridx = std::min(size, 3) - 1;
207  if (coridx == 2) {
208  if (pt > 4.5) {
209  coridx += 1;
210  }
211  if (pt > 18.) {
212  coridx += 1;
213  }
214  }
215  if (!iseb) {
216  coridx += 5;
217  }
218 
219  const GBRForestD &meanforest = *forestH_mean[coridx].product();
220  const GBRForestD &sigmaforest = *forestH_sigma[coridx].product();
221 
222  //fill array for forest evaluation
223  eval[0] = evale;
224  eval[1] = ietaix;
225  eval[2] = iphiiy;
226  if (!iseb) {
227  eval[3] = ePS1 * invE;
228  eval[4] = ePS2 * invE;
229  }
230 
231  //these are the actual BDT responses
232  double rawmean = meanforest.GetResponse(eval.data());
233  double rawsigma = sigmaforest.GetResponse(eval.data());
234 
235  //apply transformation to limited output range (matching the training)
236  double mean = iseb ? meanoffsetEB_ + meanscaleEB_ * vdt::fast_sin(rawmean)
237  : meanoffsetEE_ + meanscaleEE_ * vdt::fast_sin(rawmean);
238  double sigma = iseb ? sigmaoffsetEB_ + sigmascaleEB_ * vdt::fast_sin(rawsigma)
239  : sigmaoffsetEE_ + sigmascaleEE_ * vdt::fast_sin(rawsigma);
240 
241  //regression target is ln(Etrue/Eraw)
242  //so corrected energy is ecor=exp(mean)*e, uncertainty is exp(mean)*eraw*sigma=ecor*sigma
243  double ecor = vdt::fast_exp(mean) * e;
244  double sigmacor = sigma * ecor;
245 
246  cluster.setCorrectedEnergy(ecor);
248  cluster.setCorrectedEnergyUncertainty(sigmacor);
249  else
250  cluster.setCorrectedEnergyUncertainty(0.);
251  }
252  return;
253  }
254 
255  // Selective Readout Flags
258  evt.getByToken(ebSrFlagToken_, ebSrFlags);
259  evt.getByToken(eeSrFlagToken_, eeSrFlags);
260  if (!ebSrFlags.isValid() || !eeSrFlags.isValid())
261  throw cms::Exception("PFClusterEMEnergyCorrector")
262  << "This version of PFCluster corrections requires the SrFlagCollection information to proceed!\n";
263 
264  const unsigned int ncor = condnames_mean_.size();
265 
266  std::vector<edm::ESHandle<GBRForestD> > forestH_mean(ncor);
267  std::vector<edm::ESHandle<GBRForestD> > forestH_sigma(ncor);
268 
269  for (unsigned int icor = 0; icor < ncor; ++icor) {
270  es.get<GBRDWrapperRcd>().get(condnames_mean_[icor], forestH_mean[icor]);
271  es.get<GBRDWrapperRcd>().get(condnames_sigma_[icor], forestH_sigma[icor]);
272  }
273 
274  std::array<float, 6> evalEB;
275  std::array<float, 5> evalEE;
276 
277  for (unsigned int idx = 0; idx < cs.size(); ++idx) {
278  reco::PFCluster &cluster = cs[idx];
279  bool iseb = cluster.layer() == PFLayer::ECAL_BARREL;
280  float ePS1 = 0., ePS2 = 0.;
281  if (!iseb)
282  getAssociatedPSEnergy(idx, assoc, ePS1, ePS2);
283 
284  double e = cluster.energy();
285  double pt = cluster.pt();
286  double evale = e;
288  evale *= maxPtForMVAEvaluation_ / pt;
289  }
290  double invE = (e == 0.) ? 0. : 1. / e; //guard against dividing by 0.
291  int size = lazyTool.n5x5(cluster);
292  int reducedHits = size;
293  if (size >= 3)
294  reducedHits = 3;
295 
296  int ietaix = 0;
297  int iphiiy = 0;
298  if (iseb) {
299  EBDetId ebseed(cluster.seed());
300  ietaix = ebseed.ieta();
301  iphiiy = ebseed.iphi();
302  } else {
303  EEDetId eeseed(cluster.seed());
304  ietaix = eeseed.ix();
305  iphiiy = eeseed.iy();
306  }
307 
308  // Hardcoded number are positions of modules boundaries of ECAL
309  int signeta = (ietaix > 0) ? 1 : -1;
310  int ietamod20 = (std::abs(ietaix) < 26) ? ietaix - signeta : (ietaix - 26 * signeta) % 20;
311  int iphimod20 = (iphiiy - 1) % 20;
312 
313  int clusFlag = 0;
314  if (iseb) {
315  auto ecalUnit = readoutTool.readOutUnitOf(static_cast<EBDetId>(cluster.seed()));
316  EBSrFlagCollection::const_iterator srf = ebSrFlags->find(ecalUnit);
317  if (srf != ebSrFlags->end())
318  clusFlag = srf->value();
319  else
320  clusFlag = 3;
321  } else {
322  auto ecalUnit = readoutTool.readOutUnitOf(static_cast<EEDetId>(cluster.seed()));
323  EESrFlagCollection::const_iterator srf = eeSrFlags->find(ecalUnit);
324  if (srf != eeSrFlags->end())
325  clusFlag = srf->value();
326  else
327  clusFlag = 3;
328  }
329 
330  //find index of corrections (0-3 for EB, 4-7 for EE, depending on cluster size and raw pt)
331  int coridx = 0;
332  int regind = 0;
333  if (!iseb)
334  regind = 4;
335 
345  int ZS_bit = clusFlag >> 0 & 1;
346  int FR_bit = clusFlag >> 1 & 1;
347 
348  if (ZS_bit != 0 && FR_bit == 0)
349  coridx = 0 + regind;
350  else {
351  if (pt < 2.5)
352  coridx = 1 + regind;
353  else if (pt >= 2.5 && pt < 6.)
354  coridx = 2 + regind;
355  else if (pt >= 6.)
356  coridx = 3 + regind;
357  }
358  if (ZS_bit == 0 || clusFlag > 7) {
359  edm::LogWarning("PFClusterEMEnergyCorrector")
360  << "We can only correct regions readout in ZS (flag 1,5) or FULL readout (flag 3,7). Flag " << clusFlag
361  << " is not recognized."
362  << "\n"
363  << "Assuming FULL readout and continuing";
364  }
365 
366  const GBRForestD &meanforest = *forestH_mean[coridx].product();
367  const GBRForestD &sigmaforest = *forestH_sigma[coridx].product();
368 
369  //fill array for forest evaluation
370  if (iseb) {
371  evalEB[0] = evale;
372  evalEB[1] = ietaix;
373  evalEB[2] = iphiiy;
374  evalEB[3] = ietamod20;
375  evalEB[4] = iphimod20;
376  evalEB[5] = reducedHits;
377  } else {
378  evalEE[0] = evale;
379  evalEE[1] = ietaix;
380  evalEE[2] = iphiiy;
381  evalEE[3] = (ePS1 + ePS2) * invE;
382  evalEE[4] = reducedHits;
383  }
384 
385  //these are the actual BDT responses
386  double rawmean = 1;
387  double rawsigma = 0;
388 
389  if (iseb) {
390  rawmean = meanforest.GetResponse(evalEB.data());
391  rawsigma = sigmaforest.GetResponse(evalEB.data());
392  } else {
393  rawmean = meanforest.GetResponse(evalEE.data());
394  rawsigma = sigmaforest.GetResponse(evalEE.data());
395  }
396 
397  //apply transformation to limited output range (matching the training)
398  //the training was done with different transformations for EB and EE (width only)
399  //makes a the code a bit more cumbersome, but it is not a problem per se
400  double mean = iseb ? meanoffsetEB_ + meanscaleEB_ * vdt::fast_sin(rawmean)
401  : meanoffsetEE_ + meanscaleEE_ * vdt::fast_sin(rawmean);
402  double sigma = iseb ? sigmaoffsetEB_ + sigmascaleEB_ * vdt::fast_sin(rawsigma)
403  : sigmaoffsetEE_ + sigmascaleEE_ * vdt::fast_sin(rawsigma);
404 
405  //regression target is ln(Etrue/Eraw)
406  //so corrected energy is ecor=exp(mean)*e, uncertainty is exp(mean)*eraw*sigma=ecor*sigma
407  double ecor = iseb ? vdt::fast_exp(mean) * e : vdt::fast_exp(mean) * (e + ePS1 + ePS2);
408  double sigmacor = sigma * ecor;
409 
410  LogDebug("PFClusterEMEnergyCorrector")
411  << "ieta : iphi : ietamod20 : iphimod20 : size : reducedHits = " << ietaix << " " << iphiiy << " " << ietamod20
412  << " " << iphimod20 << " " << size << " " << reducedHits << "\n"
413  << "isEB : eraw : ePS1 : ePS2 : (eps1+eps2)/raw : Flag = " << iseb << " " << evale << " " << ePS1 << " " << ePS2
414  << " " << (ePS1 + ePS2) / evale << " " << clusFlag << "\n"
415  << "response : correction = " << exp(mean) << " " << ecor;
416 
417  cluster.setCorrectedEnergy(ecor);
419  cluster.setCorrectedEnergyUncertainty(sigmacor);
420  else
421  cluster.setCorrectedEnergyUncertainty(0.);
422  }
423 }
#define LogDebug(id)
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:100
size
Write out results.
double GetResponse(const float *vector) const
Definition: GBRForestD.h:52
int ix() const
Definition: EEDetId.h:77
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:46
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
unique_ptr< ClusterSequence > cs
std::vector< T >::const_iterator const_iterator
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:129
edm::EDGetTokenT< EESrFlagCollection > eeSrFlagToken_
std::unique_ptr< PFEnergyCalibration > calibrator_
edm::EDGetTokenT< unsigned int > bunchSpacing_
edm::EDGetTokenT< EcalRecHitCollection > recHitsEE_
std::vector< std::string > condnames_mean_50ns_
std::vector< std::string > condnames_sigma_25ns_
std::vector< std::string > condnames_mean_25ns_
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:136
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< std::string > condnames_sigma_
bool isValid() const
Definition: HandleBase.h:70
double energy() const
cluster energy
Definition: PFCluster.h:78
const_iterator end() const
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:218
edm::EDGetTokenT< EcalRecHitCollection > recHitsEB_
iterator find(key_type k)
T get() const
Definition: EventSetup.h:73
edm::EDGetTokenT< EBSrFlagCollection > ebSrFlagToken_
std::vector< std::string > condnames_sigma_50ns_
void setCorrectedEnergyUncertainty(float energyerr)
Definition: CaloCluster.h:137
void getAssociatedPSEnergy(const size_t clusIdx, const reco::PFCluster::EEtoPSAssociation &assoc, float &e1, float &e2)
std::vector< std::string > condnames_mean_
void PFClusterEMEnergyCorrector::getAssociatedPSEnergy ( const size_t  clusIdx,
const reco::PFCluster::EEtoPSAssociation assoc,
float &  e1,
float &  e2 
)
private

Definition at line 425 of file PFClusterEMEnergyCorrector.cc.

References reco::PFCluster::energy(), reco::PFCluster::layer(), PFLayer::PS1, PFLayer::PS2, and sortByKey().

Referenced by correctEnergies().

428  {
429  e1 = 0;
430  e2 = 0;
431  auto ee_key_val = std::make_pair(clusIdx, edm::Ptr<reco::PFCluster>());
432  const auto clustops = std::equal_range(assoc.begin(), assoc.end(), ee_key_val, sortByKey);
433  for (auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
434  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
435  switch (psclus->layer()) {
436  case PFLayer::PS1:
437  e1 += psclus->energy();
438  break;
439  case PFLayer::PS2:
440  e2 += psclus->energy();
441  break;
442  default:
443  break;
444  }
445  }
446 }
bool sortByKey(const EEPSPair &a, const EEPSPair &b)
PFClusterEMEnergyCorrector& PFClusterEMEnergyCorrector::operator= ( const PFClusterEMEnergyCorrector )
delete

Member Data Documentation

bool PFClusterEMEnergyCorrector::applyCrackCorrections_
private

Definition at line 57 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::applyMVACorrections_
private

Definition at line 58 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::autoDetectBunchSpacing_
private

Definition at line 61 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

edm::EDGetTokenT<unsigned int> PFClusterEMEnergyCorrector::bunchSpacing_
private

Definition at line 46 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

int PFClusterEMEnergyCorrector::bunchSpacingManual_
private

Definition at line 62 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

std::unique_ptr<PFEnergyCalibration> PFClusterEMEnergyCorrector::calibrator_
private

Definition at line 64 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies().

std::vector<std::string> PFClusterEMEnergyCorrector::condnames_mean_
private

Definition at line 48 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

std::vector<std::string> PFClusterEMEnergyCorrector::condnames_mean_25ns_
private

Definition at line 51 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

std::vector<std::string> PFClusterEMEnergyCorrector::condnames_mean_50ns_
private

Definition at line 53 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

std::vector<std::string> PFClusterEMEnergyCorrector::condnames_sigma_
private

Definition at line 49 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

std::vector<std::string> PFClusterEMEnergyCorrector::condnames_sigma_25ns_
private

Definition at line 52 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

std::vector<std::string> PFClusterEMEnergyCorrector::condnames_sigma_50ns_
private

Definition at line 54 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

edm::EDGetTokenT<EBSrFlagCollection> PFClusterEMEnergyCorrector::ebSrFlagToken_
private

Definition at line 40 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

edm::EDGetTokenT<EESrFlagCollection> PFClusterEMEnergyCorrector::eeSrFlagToken_
private

Definition at line 41 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::maxPtForMVAEvaluation_
private

Definition at line 38 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimhighEB_
private

Definition at line 71 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimhighEE_
private

Definition at line 76 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimlowEB_
private

Definition at line 70 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimlowEE_
private

Definition at line 75 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanoffsetEB_
private

Definition at line 72 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanoffsetEE_
private

Definition at line 77 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanscaleEB_
private

Definition at line 73 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanscaleEE_
private

Definition at line 78 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

edm::EDGetTokenT<EcalRecHitCollection> PFClusterEMEnergyCorrector::recHitsEB_
private

Definition at line 44 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

edm::EDGetTokenT<EcalRecHitCollection> PFClusterEMEnergyCorrector::recHitsEE_
private

Definition at line 45 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::setEnergyUncertainty_
private

Definition at line 59 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimhighEB_
private

Definition at line 81 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimhighEE_
private

Definition at line 86 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimlowEB_
private

Definition at line 80 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimlowEE_
private

Definition at line 85 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmaoffsetEB_
private

Definition at line 82 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmaoffsetEE_
private

Definition at line 87 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmascaleEB_
private

Definition at line 83 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmascaleEE_
private

Definition at line 88 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::srfAwareCorrection_
private

Definition at line 56 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().