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 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 Attributes

bool _applyCrackCorrections
 
bool _applyMVACorrections
 
std::unique_ptr
< PFEnergyCalibration
_calibrator
 
std::vector< std::string > _condnames_mean_25ns
 
std::vector< std::string > _condnames_mean_50ns
 
std::vector< std::string > _condnames_sigma_25ns
 
std::vector< std::string > _condnames_sigma_50ns
 
double _maxPtForMVAEvaluation
 
edm::EDGetTokenT
< EcalRecHitCollection
_recHitsEB
 
edm::EDGetTokenT
< EcalRecHitCollection
_recHitsEE
 
bool autoDetectBunchSpacing_
 
edm::EDGetTokenT< int > bunchSpacing_
 
int bunchSpacingManual_
 

Detailed Description

Definition at line 17 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, _condnames_mean_25ns, _condnames_mean_50ns, _condnames_sigma_25ns, _condnames_sigma_50ns, _maxPtForMVAEvaluation, _recHitsEB, _recHitsEE, autoDetectBunchSpacing_, bunchSpacing_, bunchSpacingManual_, edm::ParameterSet::getParameter(), and HLT_25ns14e33_v1_cff::InputTag.

16  :
18 
19  _applyCrackCorrections = conf.getParameter<bool>("applyCrackCorrections");
20  _applyMVACorrections = conf.getParameter<bool>("applyMVACorrections");
21  _maxPtForMVAEvaluation = conf.getParameter<double>("maxPtForMVAEvaluation");
22 
23 
24  if (_applyMVACorrections) {
25  _recHitsEB = cc.consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitsEBLabel"));
26  _recHitsEE = cc.consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitsEELabel"));
27 
28  autoDetectBunchSpacing_ = conf.getParameter<bool>("autoDetectBunchSpacing");
29 
31  bunchSpacing_ = cc.consumes<int>(edm::InputTag("addPileupInfo","bunchSpacing"));
33  }
34  else {
35  bunchSpacingManual_ = conf.getParameter<int>("bunchSpacing");
36  }
37 
38  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EB_pfSize1_mean_50ns");
39  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EB_pfSize2_mean_50ns");
40  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin1_mean_50ns");
41  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin2_mean_50ns");
42  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin3_mean_50ns");
43  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EE_pfSize1_mean_50ns");
44  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EE_pfSize2_mean_50ns");
45  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin1_mean_50ns");
46  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin2_mean_50ns");
47  _condnames_mean_50ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin3_mean_50ns");
48 
49  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EB_pfSize1_sigma_50ns");
50  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EB_pfSize2_sigma_50ns");
51  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin1_sigma_50ns");
52  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin2_sigma_50ns");
53  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin3_sigma_50ns");
54  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EE_pfSize1_sigma_50ns");
55  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EE_pfSize2_sigma_50ns");
56  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin1_sigma_50ns");
57  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin2_sigma_50ns");
58  _condnames_sigma_50ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin3_sigma_50ns");
59 
60  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EB_pfSize1_mean_25ns");
61  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EB_pfSize2_mean_25ns");
62  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin1_mean_25ns");
63  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin2_mean_25ns");
64  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin3_mean_25ns");
65  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EE_pfSize1_mean_25ns");
66  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EE_pfSize2_mean_25ns");
67  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin1_mean_25ns");
68  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin2_mean_25ns");
69  _condnames_mean_25ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin3_mean_25ns");
70 
71  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EB_pfSize1_sigma_25ns");
72  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EB_pfSize2_sigma_25ns");
73  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin1_sigma_25ns");
74  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin2_sigma_25ns");
75  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EB_pfSize3_ptbin3_sigma_25ns");
76  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EE_pfSize1_sigma_25ns");
77  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EE_pfSize2_sigma_25ns");
78  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin1_sigma_25ns");
79  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin2_sigma_25ns");
80  _condnames_sigma_25ns.push_back("ecalPFClusterCorV2_EE_pfSize3_ptbin3_sigma_25ns");
81  }
82 
83 
84 
85 }
std::vector< std::string > _condnames_mean_50ns
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
std::vector< std::string > _condnames_sigma_25ns
std::vector< std::string > _condnames_mean_25ns
edm::EDGetTokenT< int > bunchSpacing_
edm::EDGetTokenT< EcalRecHitCollection > _recHitsEB
std::unique_ptr< PFEnergyCalibration > _calibrator
edm::EDGetTokenT< EcalRecHitCollection > _recHitsEE
std::vector< std::string > _condnames_sigma_50ns
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 
)

Definition at line 88 of file PFClusterEMEnergyCorrector.cc.

References _applyCrackCorrections, _applyMVACorrections, _calibrator, _condnames_mean_25ns, _condnames_mean_50ns, _condnames_sigma_25ns, _condnames_sigma_50ns, _maxPtForMVAEvaluation, _recHitsEB, _recHitsEE, autoDetectBunchSpacing_, bunchSpacing_, bunchSpacingManual_, alignCSCRings::e, PFLayer::ECAL_BARREL, reco::PFCluster::energy(), create_public_lumi_plots::exp, edm::EventSetup::get(), edm::Event::getByToken(), GBRForestD::GetResponse(), customizeTrackingMonitorSeedNumber::idx, EBDetId::ieta(), edm::EventBase::isRealData(), EEDetId::ix(), reco::PFCluster::layer(), timingPdfMaker::mean, min(), PFLayer::PS1, PFLayer::PS2, EnergyCorrector::pt, reco::PFCluster::pt(), DTTTrigCorrFirst::run, edm::Event::run(), reco::CaloCluster::seed(), reco::CaloCluster::setCorrectedEnergy(), reco::CaloCluster::setCorrectedEnergyUncertainty(), and findQualityFiles::size.

88  {
89 
90  //legacy corrections
91  if (!_applyMVACorrections) {
92  for (unsigned int idx = 0; idx<cs.size(); ++idx) {
93  reco::PFCluster &cluster = cs[idx];
94  bool iseb = cluster.layer() == PFLayer::ECAL_BARREL;
95 
96  //compute preshower energies for endcap clusters
97  double ePS1=0, ePS2=0;
98  if(!iseb) {
99  auto ee_key_val = std::make_pair(idx,edm::Ptr<reco::PFCluster>());
100  const auto clustops = std::equal_range(assoc.begin(),
101  assoc.end(),
102  ee_key_val,
103  sortByKey);
104  for( auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
105  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
106  switch( psclus->layer() ) {
107  case PFLayer::PS1:
108  ePS1 += psclus->energy();
109  break;
110  case PFLayer::PS2:
111  ePS2 += psclus->energy();
112  break;
113  default:
114  break;
115  }
116  }
117  }
118 
119  double correctedEnergy = _calibrator->energyEm(cluster,ePS1,ePS2,_applyCrackCorrections);
120  cluster.setCorrectedEnergy(correctedEnergy);
121 
122  }
123  return;
124  }
125 
126  int bunchspacing = 450;
127 
129  if (evt.isRealData()) {
130  edm::RunNumber_t run = evt.run();
131  if (run == 178003 ||
132  run == 178004 ||
133  run == 209089 ||
134  run == 209106 ||
135  run == 209109 ||
136  run == 209146 ||
137  run == 209148 ||
138  run == 209151) {
139  bunchspacing = 25;
140  }
141  else if (run < 253000) {
142  bunchspacing = 50;
143  }
144  else {
145  bunchspacing = 25;
146  }
147  }
148  else {
149  edm::Handle<int> bunchSpacingH;
150  evt.getByToken(bunchSpacing_,bunchSpacingH);
151  bunchspacing = *bunchSpacingH;
152  }
153  }
154  else {
155  bunchspacing = bunchSpacingManual_;
156  }
157 
158  const std::vector<std::string> condnames_mean = (bunchspacing == 25) ? _condnames_mean_25ns : _condnames_mean_50ns;
159  const std::vector<std::string> condnames_sigma = (bunchspacing == 25) ? _condnames_sigma_25ns : _condnames_sigma_50ns;
160 
161  const unsigned int ncor = condnames_mean.size();
162 
163  std::vector<edm::ESHandle<GBRForestD> > forestH_mean(ncor);
164  std::vector<edm::ESHandle<GBRForestD> > forestH_sigma(ncor);
165 
166  for (unsigned int icor=0; icor<ncor; ++icor) {
167  es.get<GBRDWrapperRcd>().get(condnames_mean[icor],forestH_mean[icor]);
168  es.get<GBRDWrapperRcd>().get(condnames_sigma[icor],forestH_sigma[icor]);
169  }
170 
171  std::array<float,5> eval;
172 
173  EcalClusterLazyTools lazyTool(evt, es, _recHitsEB, _recHitsEE);
174 
175  //magic numbers for MINUIT-like transformation of BDT output onto limited range
176  //(These should be stored inside the conditions object in the future as well)
177  const double meanlimlow = -0.336;
178  const double meanlimhigh = 0.916;
179  const double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
180  const double meanscale = 0.5*(meanlimhigh-meanlimlow);
181 
182  const double sigmalimlow = 0.001;
183  const double sigmalimhigh = 0.4;
184  const double sigmaoffset = sigmalimlow + 0.5*(sigmalimhigh-sigmalimlow);
185  const double sigmascale = 0.5*(sigmalimhigh-sigmalimlow);
186 
187  for (unsigned int idx = 0; idx<cs.size(); ++idx) {
188 
189  reco::PFCluster &cluster = cs[idx];
190 
191  double e = cluster.energy();
192  double pt = cluster.pt();
193 
194  //limit raw energy value used to evaluate corrections
195  //to avoid bad extrapolation
196  double evale = e;
198  evale *= _maxPtForMVAEvaluation/pt;
199  }
200 
201  double invE = (e == 0.) ? 0. : 1./e; //guard against dividing by 0.
202 
203  int size = lazyTool.n5x5(cluster);
204 
205  bool iseb = cluster.layer() == PFLayer::ECAL_BARREL;
206 
207  //find index of corrections (0-4 for EB, 5-9 for EE, depending on cluster size and raw pt)
208  int coridx = std::min(size,3)-1;
209  if (coridx==2) {
210  if (pt>4.5) {
211  coridx += 1;
212  }
213  if (pt>18.) {
214  coridx += 1;
215  }
216  }
217  if (!iseb) {
218  coridx += 5;
219  }
220 
221  const GBRForestD &meanforest = *forestH_mean[coridx].product();
222  const GBRForestD &sigmaforest = *forestH_sigma[coridx].product();
223 
224  //find seed crystal indices
225  int ietaix=0;
226  int iphiiy=0;
227  if (iseb) {
228  EBDetId ebseed(cluster.seed());
229  ietaix = ebseed.ieta();
230  iphiiy = ebseed.iphi();
231  }
232  else {
233  EEDetId eeseed(cluster.seed());
234  ietaix = eeseed.ix();
235  iphiiy = eeseed.iy();
236  }
237 
238 
239  //compute preshower energies for endcap clusters
240  double ePS1=0, ePS2=0;
241  if(!iseb) {
242  auto ee_key_val = std::make_pair(idx,edm::Ptr<reco::PFCluster>());
243  const auto clustops = std::equal_range(assoc.begin(),
244  assoc.end(),
245  ee_key_val,
246  sortByKey);
247  for( auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
248  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
249  switch( psclus->layer() ) {
250  case PFLayer::PS1:
251  ePS1 += psclus->energy();
252  break;
253  case PFLayer::PS2:
254  ePS2 += psclus->energy();
255  break;
256  default:
257  break;
258  }
259  }
260  }
261 
262  //fill array for forest evaluation
263  eval[0] = evale;
264  eval[1] = ietaix;
265  eval[2] = iphiiy;
266  if (!iseb) {
267  eval[3] = ePS1*invE;
268  eval[4] = ePS2*invE;
269  }
270 
271  //these are the actual BDT responses
272  double rawmean = meanforest.GetResponse(eval.data());
273  double rawsigma = sigmaforest.GetResponse(eval.data());
274 
275  //apply transformation to limited output range (matching the training)
276  double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
277  double sigma = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma);
278 
279  //regression target is ln(Etrue/Eraw)
280  //so corrected energy is ecor=exp(mean)*e, uncertainty is exp(mean)*eraw*sigma=ecor*sigma
281  double ecor = exp(mean)*e;
282  double sigmacor = sigma*ecor;
283 
284  cluster.setCorrectedEnergy(ecor);
285  cluster.setCorrectedEnergyUncertainty(sigmacor);
286 
287  }
288 
289 }
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:86
std::vector< std::string > _condnames_mean_50ns
double GetResponse(const float *vector) const
Definition: GBRForestD.h:51
auto_ptr< ClusterSequence > cs
int ix() const
Definition: EEDetId.h:76
std::vector< std::string > _condnames_sigma_25ns
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:138
std::vector< std::string > _condnames_mean_25ns
bool isRealData() const
Definition: EventBase.h:64
edm::EDGetTokenT< int > bunchSpacing_
edm::EDGetTokenT< EcalRecHitCollection > _recHitsEB
std::unique_ptr< PFEnergyCalibration > _calibrator
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:109
RunNumber_t run() const
Definition: Event.h:87
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
T min(T a, T b)
Definition: MathUtil.h:58
double energy() const
cluster energy
Definition: PFCluster.h:82
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:202
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
const T & get() const
Definition: EventSetup.h:55
edm::EDGetTokenT< EcalRecHitCollection > _recHitsEE
std::vector< std::string > _condnames_sigma_50ns
unsigned int RunNumber_t
void setCorrectedEnergyUncertainty(float energyerr)
Definition: CaloCluster.h:110
tuple size
Write out results.
PFClusterEMEnergyCorrector& PFClusterEMEnergyCorrector::operator= ( const PFClusterEMEnergyCorrector )
delete

Member Data Documentation

bool PFClusterEMEnergyCorrector::_applyCrackCorrections
private

Definition at line 26 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::_applyMVACorrections
private

Definition at line 27 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 43 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies().

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

Definition at line 40 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 38 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 41 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 39 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

double PFClusterEMEnergyCorrector::_maxPtForMVAEvaluation
private

Definition at line 28 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 35 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 36 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::autoDetectBunchSpacing_
private

Definition at line 30 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 33 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

int PFClusterEMEnergyCorrector::bunchSpacingManual_
private

Definition at line 31 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().