CMS 3D CMS Logo

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< unsigned 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_, and edm::ParameterSet::getParameter().

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<unsigned int>(edm::InputTag("bunchSpacingProducer"));
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< unsigned 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_, MillePedeFileConverter_cfg::e, PFLayer::ECAL_BARREL, reco::PFCluster::energy(), JetChargeProducer_cfi::exp, edm::EventSetup::get(), edm::Event::getByToken(), GBRForestD::GetResponse(), training_settings::idx, EBDetId::ieta(), EEDetId::ix(), reco::PFCluster::layer(), RecoTauDiscriminantConfiguration::mean, min(), PFLayer::PS1, PFLayer::PS2, EnergyCorrector::pt, reco::PFCluster::pt(), reco::CaloCluster::seed(), reco::CaloCluster::setCorrectedEnergy(), reco::CaloCluster::setCorrectedEnergyUncertainty(), ctppsDiamondLocalTracks_cfi::sigma, findQualityFiles::size, and sortByKey().

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  edm::Handle<unsigned int> bunchSpacingH;
130  evt.getByToken(bunchSpacing_,bunchSpacingH);
131  bunchspacing = *bunchSpacingH;
132  }
133  else {
134  bunchspacing = bunchSpacingManual_;
135  }
136 
137  const std::vector<std::string> condnames_mean = (bunchspacing == 25) ? _condnames_mean_25ns : _condnames_mean_50ns;
138  const std::vector<std::string> condnames_sigma = (bunchspacing == 25) ? _condnames_sigma_25ns : _condnames_sigma_50ns;
139 
140  const unsigned int ncor = condnames_mean.size();
141 
142  std::vector<edm::ESHandle<GBRForestD> > forestH_mean(ncor);
143  std::vector<edm::ESHandle<GBRForestD> > forestH_sigma(ncor);
144 
145  for (unsigned int icor=0; icor<ncor; ++icor) {
146  es.get<GBRDWrapperRcd>().get(condnames_mean[icor],forestH_mean[icor]);
147  es.get<GBRDWrapperRcd>().get(condnames_sigma[icor],forestH_sigma[icor]);
148  }
149 
150  std::array<float,5> eval;
151 
152  EcalClusterLazyTools lazyTool(evt, es, _recHitsEB, _recHitsEE);
153 
154  //magic numbers for MINUIT-like transformation of BDT output onto limited range
155  //(These should be stored inside the conditions object in the future as well)
156  const double meanlimlow = -0.336;
157  const double meanlimhigh = 0.916;
158  const double meanoffset = meanlimlow + 0.5*(meanlimhigh-meanlimlow);
159  const double meanscale = 0.5*(meanlimhigh-meanlimlow);
160 
161  const double sigmalimlow = 0.001;
162  const double sigmalimhigh = 0.4;
163  const double sigmaoffset = sigmalimlow + 0.5*(sigmalimhigh-sigmalimlow);
164  const double sigmascale = 0.5*(sigmalimhigh-sigmalimlow);
165 
166  for (unsigned int idx = 0; idx<cs.size(); ++idx) {
167 
168  reco::PFCluster &cluster = cs[idx];
169 
170  double e = cluster.energy();
171  double pt = cluster.pt();
172 
173  //limit raw energy value used to evaluate corrections
174  //to avoid bad extrapolation
175  double evale = e;
177  evale *= _maxPtForMVAEvaluation/pt;
178  }
179 
180  double invE = (e == 0.) ? 0. : 1./e; //guard against dividing by 0.
181 
182  int size = lazyTool.n5x5(cluster);
183 
184  bool iseb = cluster.layer() == PFLayer::ECAL_BARREL;
185 
186  //find index of corrections (0-4 for EB, 5-9 for EE, depending on cluster size and raw pt)
187  int coridx = std::min(size,3)-1;
188  if (coridx==2) {
189  if (pt>4.5) {
190  coridx += 1;
191  }
192  if (pt>18.) {
193  coridx += 1;
194  }
195  }
196  if (!iseb) {
197  coridx += 5;
198  }
199 
200  const GBRForestD &meanforest = *forestH_mean[coridx].product();
201  const GBRForestD &sigmaforest = *forestH_sigma[coridx].product();
202 
203  //find seed crystal indices
204  int ietaix=0;
205  int iphiiy=0;
206  if (iseb) {
207  EBDetId ebseed(cluster.seed());
208  ietaix = ebseed.ieta();
209  iphiiy = ebseed.iphi();
210  }
211  else {
212  EEDetId eeseed(cluster.seed());
213  ietaix = eeseed.ix();
214  iphiiy = eeseed.iy();
215  }
216 
217 
218  //compute preshower energies for endcap clusters
219  double ePS1=0, ePS2=0;
220  if(!iseb) {
221  auto ee_key_val = std::make_pair(idx,edm::Ptr<reco::PFCluster>());
222  const auto clustops = std::equal_range(assoc.begin(),
223  assoc.end(),
224  ee_key_val,
225  sortByKey);
226  for( auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
227  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
228  switch( psclus->layer() ) {
229  case PFLayer::PS1:
230  ePS1 += psclus->energy();
231  break;
232  case PFLayer::PS2:
233  ePS2 += psclus->energy();
234  break;
235  default:
236  break;
237  }
238  }
239  }
240 
241  //fill array for forest evaluation
242  eval[0] = evale;
243  eval[1] = ietaix;
244  eval[2] = iphiiy;
245  if (!iseb) {
246  eval[3] = ePS1*invE;
247  eval[4] = ePS2*invE;
248  }
249 
250  //these are the actual BDT responses
251  double rawmean = meanforest.GetResponse(eval.data());
252  double rawsigma = sigmaforest.GetResponse(eval.data());
253 
254  //apply transformation to limited output range (matching the training)
255  double mean = meanoffset + meanscale*vdt::fast_sin(rawmean);
256  double sigma = sigmaoffset + sigmascale*vdt::fast_sin(rawsigma);
257 
258  //regression target is ln(Etrue/Eraw)
259  //so corrected energy is ecor=exp(mean)*e, uncertainty is exp(mean)*eraw*sigma=ecor*sigma
260  double ecor = exp(mean)*e;
261  double sigmacor = sigma*ecor;
262 
263  cluster.setCorrectedEnergy(ecor);
264  cluster.setCorrectedEnergyUncertainty(sigmacor);
265 
266  }
267 
268 }
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:128
std::vector< std::string > _condnames_mean_50ns
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
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:460
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:139
std::vector< std::string > _condnames_mean_25ns
edm::EDGetTokenT< unsigned int > bunchSpacing_
edm::EDGetTokenT< EcalRecHitCollection > _recHitsEB
std::unique_ptr< PFEnergyCalibration > _calibrator
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:112
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:205
const T & get() const
Definition: EventSetup.h:56
bool sortByKey(const EEPSPair &a, const EEPSPair &b)
edm::EDGetTokenT< EcalRecHitCollection > _recHitsEE
std::vector< std::string > _condnames_sigma_50ns
void setCorrectedEnergyUncertainty(float energyerr)
Definition: CaloCluster.h:113
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<unsigned 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().