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_
 
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 16 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(), maxPtForMVAEvaluation_, meanlimhighEB_, meanlimhighEE_, meanlimlowEB_, meanlimlowEE_, meanoffsetEB_, meanoffsetEE_, meanscaleEB_, meanscaleEE_, recHitsEB_, recHitsEE_, sigmalimhighEB_, sigmalimhighEE_, sigmalimlowEB_, sigmalimlowEE_, sigmaoffsetEB_, sigmaoffsetEE_, sigmascaleEB_, sigmascaleEE_, and srfAwareCorrection_.

16  :
18 
19  applyCrackCorrections_ = conf.getParameter<bool>("applyCrackCorrections");
20  applyMVACorrections_ = conf.getParameter<bool>("applyMVACorrections");
21  srfAwareCorrection_ = conf.getParameter<bool>("srfAwareCorrection");
22 
23  maxPtForMVAEvaluation_ = conf.getParameter<double>("maxPtForMVAEvaluation");
24 
25  if (applyMVACorrections_) {
26 
27  meanlimlowEB_ = -0.336;
28  meanlimhighEB_ = 0.916;
31 
32  meanlimlowEE_ = -0.336;
33  meanlimhighEE_ = 0.916;
36 
37  sigmalimlowEB_ = 0.001;
38  sigmalimhighEB_ = 0.4;
41 
42  sigmalimlowEE_ = 0.001;
43  sigmalimhighEE_ = 0.4;
46 
47  recHitsEB_ = cc.consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitsEBLabel"));
48  recHitsEE_ = cc.consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitsEELabel"));
49  autoDetectBunchSpacing_ = conf.getParameter<bool>("autoDetectBunchSpacing");
50 
52  bunchSpacing_ = cc.consumes<unsigned int>(edm::InputTag("bunchSpacingProducer"));
54  } else {
55  bunchSpacingManual_ = conf.getParameter<int>("bunchSpacing");
56  }
57 
59  "ecalPFClusterCorV2_EB_pfSize1_mean_25ns",
60  "ecalPFClusterCorV2_EB_pfSize2_mean_25ns",
61  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_mean_25ns",
62  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_mean_25ns",
63  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_mean_25ns",
64  "ecalPFClusterCorV2_EE_pfSize1_mean_25ns",
65  "ecalPFClusterCorV2_EE_pfSize2_mean_25ns",
66  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_mean_25ns",
67  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_mean_25ns",
68  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_mean_25ns"});
70  "ecalPFClusterCorV2_EB_pfSize1_sigma_25ns",
71  "ecalPFClusterCorV2_EB_pfSize2_sigma_25ns",
72  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_sigma_25ns",
73  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_sigma_25ns",
74  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_sigma_25ns",
75  "ecalPFClusterCorV2_EE_pfSize1_sigma_25ns",
76  "ecalPFClusterCorV2_EE_pfSize2_sigma_25ns",
77  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_sigma_25ns",
78  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_sigma_25ns",
79  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_sigma_25ns"});
81  "ecalPFClusterCorV2_EB_pfSize1_mean_50ns",
82  "ecalPFClusterCorV2_EB_pfSize2_mean_50ns",
83  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_mean_50ns",
84  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_mean_50ns",
85  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_mean_50ns",
86  "ecalPFClusterCorV2_EE_pfSize1_mean_50ns",
87  "ecalPFClusterCorV2_EE_pfSize2_mean_50ns",
88  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_mean_50ns",
89  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_mean_50ns",
90  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_mean_50ns"});
92  "ecalPFClusterCorV2_EB_pfSize1_sigma_50ns",
93  "ecalPFClusterCorV2_EB_pfSize2_sigma_50ns",
94  "ecalPFClusterCorV2_EB_pfSize3_ptbin1_sigma_50ns",
95  "ecalPFClusterCorV2_EB_pfSize3_ptbin2_sigma_50ns",
96  "ecalPFClusterCorV2_EB_pfSize3_ptbin3_sigma_50ns",
97  "ecalPFClusterCorV2_EE_pfSize1_sigma_50ns",
98  "ecalPFClusterCorV2_EE_pfSize2_sigma_50ns",
99  "ecalPFClusterCorV2_EE_pfSize3_ptbin1_sigma_50ns",
100  "ecalPFClusterCorV2_EE_pfSize3_ptbin2_sigma_50ns",
101  "ecalPFClusterCorV2_EE_pfSize3_ptbin3_sigma_50ns"});
102 
103  if (srfAwareCorrection_) {
104 
105  sigmalimlowEE_ = 0.001;
106  sigmalimhighEE_ = 0.1;
109 
112 
113  condnames_mean_.insert(condnames_mean_.end(), {
114  "ecalPFClusterCor2017V2_EB_ZS_mean_25ns",
115  "ecalPFClusterCor2017V2_EB_Full_ptbin1_mean_25ns",
116  "ecalPFClusterCor2017V2_EB_Full_ptbin2_mean_25ns",
117  "ecalPFClusterCor2017V2_EB_Full_ptbin3_mean_25ns",
118  "ecalPFClusterCor2017V2_EE_ZS_mean_25ns",
119  "ecalPFClusterCor2017V2_EE_Full_ptbin1_mean_25ns",
120  "ecalPFClusterCor2017V2_EE_Full_ptbin2_mean_25ns",
121  "ecalPFClusterCor2017V2_EE_Full_ptbin3_mean_25ns"});
122 
123  condnames_sigma_.insert(condnames_sigma_.end(), {
124  "ecalPFClusterCor2017V2_EB_ZS_sigma_25ns",
125  "ecalPFClusterCor2017V2_EB_Full_ptbin1_sigma_25ns",
126  "ecalPFClusterCor2017V2_EB_Full_ptbin2_sigma_25ns",
127  "ecalPFClusterCor2017V2_EB_Full_ptbin3_sigma_25ns",
128  "ecalPFClusterCor2017V2_EE_ZS_sigma_25ns",
129  "ecalPFClusterCor2017V2_EE_Full_ptbin1_sigma_25ns",
130  "ecalPFClusterCor2017V2_EE_Full_ptbin2_sigma_25ns",
131  "ecalPFClusterCor2017V2_EE_Full_ptbin3_sigma_25ns"});
132  }
133  }
134 
135 }
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 137 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(), EnergyCorrector::pt, reco::PFCluster::pt(), EcalReadoutTools::readOutUnitOf(), recHitsEB_, recHitsEE_, reco::CaloCluster::seed(), reco::CaloCluster::setCorrectedEnergy(), reco::CaloCluster::setCorrectedEnergyUncertainty(), sigmaoffsetEB_, sigmaoffsetEE_, sigmascaleEB_, sigmascaleEE_, findQualityFiles::size, and srfAwareCorrection_.

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

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

Referenced by correctEnergies().

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

Member Data Documentation

bool PFClusterEMEnergyCorrector::applyCrackCorrections_
private

Definition at line 54 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::applyMVACorrections_
private

Definition at line 55 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::autoDetectBunchSpacing_
private

Definition at line 57 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 43 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

int PFClusterEMEnergyCorrector::bunchSpacingManual_
private

Definition at line 58 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 60 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies().

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

Definition at line 45 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 48 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 50 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 46 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 49 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 51 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 37 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 38 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::maxPtForMVAEvaluation_
private

Definition at line 35 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimhighEB_
private

Definition at line 64 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimhighEE_
private

Definition at line 69 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimlowEB_
private

Definition at line 63 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanlimlowEE_
private

Definition at line 68 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanoffsetEB_
private

Definition at line 65 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanoffsetEE_
private

Definition at line 70 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanscaleEB_
private

Definition at line 66 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::meanscaleEE_
private

Definition at line 71 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 41 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 42 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimhighEB_
private

Definition at line 74 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimhighEE_
private

Definition at line 79 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimlowEB_
private

Definition at line 73 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmalimlowEE_
private

Definition at line 78 of file PFClusterEMEnergyCorrector.h.

Referenced by PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmaoffsetEB_
private

Definition at line 75 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmaoffsetEE_
private

Definition at line 80 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmascaleEB_
private

Definition at line 76 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::sigmascaleEE_
private

Definition at line 81 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::srfAwareCorrection_
private

Definition at line 53 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().