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

References _applyCrackCorrections, _applyMVACorrections, _calibrator, _condnames_mean_25ns, _condnames_mean_50ns, _condnames_sigma_25ns, _condnames_sigma_50ns, _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.

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

Referenced by correctEnergies().

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

Definition at line 39 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 37 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 40 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 38 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 34 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 35 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

bool PFClusterEMEnergyCorrector::autoDetectBunchSpacing_
private

Definition at line 29 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

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

Definition at line 32 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().

int PFClusterEMEnergyCorrector::bunchSpacingManual_
private

Definition at line 30 of file PFClusterEMEnergyCorrector.h.

Referenced by correctEnergies(), and PFClusterEMEnergyCorrector().