CMS 3D CMS Logo

LowPtGsfElectronFeatures.cc
Go to the documentation of this file.
9 #include "TVector3.h"
10 #include <algorithm>
11 
12 namespace lowptgsfeleid {
13 
14  std::vector<float> features_V1(reco::GsfElectron const& ele, float rho, float unbiased, float field_z) {
15  float eid_rho = -999.;
16  float eid_sc_eta = -999.;
17  float eid_shape_full5x5_r9 = -999.;
18  float eid_sc_etaWidth = -999.;
19  float eid_sc_phiWidth = -999.;
20  float eid_shape_full5x5_HoverE = -999.;
21  float eid_trk_nhits = -999.;
22  float eid_trk_chi2red = -999.;
23  float eid_gsf_chi2red = -999.;
24  float eid_brem_frac = -999.;
25  float eid_gsf_nhits = -999.;
26  float eid_match_SC_EoverP = -999.;
27  float eid_match_eclu_EoverP = -999.;
28  float eid_match_SC_dEta = -999.;
29  float eid_match_SC_dPhi = -999.;
30  float eid_match_seed_dEta = -999.;
31  float eid_sc_E = -999.;
32  float eid_trk_p = -999.;
33  float gsf_mode_p = -999.;
34  float core_shFracHits = -999.;
35  float gsf_bdtout1 = -999.;
36  float gsf_dr = -999.;
37  float trk_dr = -999.;
38  float sc_Nclus = -999.;
39  float sc_clus1_nxtal = -999.;
40  float sc_clus1_dphi = -999.;
41  float sc_clus2_dphi = -999.;
42  float sc_clus1_deta = -999.;
43  float sc_clus2_deta = -999.;
44  float sc_clus1_E = -999.;
45  float sc_clus2_E = -999.;
46  float sc_clus1_E_ov_p = -999.;
47  float sc_clus2_E_ov_p = -999.;
48 
49  // KF tracks
50  if (ele.core().isNonnull()) {
52  if (trk.isNonnull()) {
53  eid_trk_p = (float)trk->p();
54  eid_trk_nhits = (float)trk->found();
55  eid_trk_chi2red = (float)trk->normalizedChi2();
56  TVector3 trkTV3(0, 0, 0);
57  trkTV3.SetPtEtaPhi(trk->pt(), trk->eta(), trk->phi());
58  TVector3 eleTV3(0, 0, 0);
59  eleTV3.SetPtEtaPhi(ele.pt(), ele.eta(), ele.phi());
60  trk_dr = eleTV3.DeltaR(trkTV3);
61  }
62  }
63 
64  // GSF tracks
65  if (ele.core().isNonnull()) {
67  if (gsf.isNonnull()) {
68  gsf_mode_p = gsf->pMode();
69  eid_gsf_nhits = (float)gsf->found();
70  eid_gsf_chi2red = gsf->normalizedChi2();
71  TVector3 gsfTV3(0, 0, 0);
72  gsfTV3.SetPtEtaPhi(gsf->ptMode(), gsf->etaMode(), gsf->phiMode());
73  TVector3 eleTV3(0, 0, 0);
74  eleTV3.SetPtEtaPhi(ele.pt(), ele.eta(), ele.phi());
75  gsf_dr = eleTV3.DeltaR(gsfTV3);
76  }
77  }
78 
79  // Super clusters
80  if (ele.core().isNonnull()) {
82  if (sc.isNonnull()) {
83  eid_sc_E = sc->energy();
84  eid_sc_eta = sc->eta();
85  eid_sc_etaWidth = sc->etaWidth();
86  eid_sc_phiWidth = sc->phiWidth();
87  sc_Nclus = sc->clustersSize();
88  }
89  }
90 
91  // Track-cluster matching
92  eid_match_seed_dEta = ele.deltaEtaSeedClusterTrackAtCalo();
93  eid_match_eclu_EoverP = (1. / ele.ecalEnergy()) - (1. / ele.p());
94  eid_match_SC_EoverP = ele.eSuperClusterOverP();
95  eid_match_SC_dEta = ele.deltaEtaSuperClusterTrackAtVtx();
96  eid_match_SC_dPhi = ele.deltaPhiSuperClusterTrackAtVtx();
97 
98  // Shower shape vars
99  eid_shape_full5x5_HoverE = ele.full5x5_hcalOverEcal();
100  eid_shape_full5x5_r9 = ele.full5x5_r9();
101 
102  // Misc
103  eid_rho = rho;
104 
105  eid_brem_frac = ele.fbrem();
106  core_shFracHits = ele.shFracInnerHits();
107 
108  // Unbiased BDT from ElectronSeed
109  gsf_bdtout1 = unbiased;
110 
111  // Clusters
112  if (ele.core().isNonnull()) {
114  if (gsf.isNonnull()) {
116  if (sc.isNonnull()) {
117  // Propagate electron track to ECAL surface
118  double mass2 = 0.000511 * 0.000511;
119  float p2 = pow(gsf->p(), 2);
120  float energy = sqrt(mass2 + p2);
121  XYZTLorentzVector mom = XYZTLorentzVector(gsf->px(), gsf->py(), gsf->pz(), energy);
122  XYZTLorentzVector pos = XYZTLorentzVector(gsf->vx(), gsf->vy(), gsf->vz(), 0.);
123  BaseParticlePropagator propagator(RawParticle(mom, pos, gsf->charge()), 0, 0, field_z);
124 
125  propagator.propagateToEcalEntrance(true); // true only first half loop , false more than one loop
126  bool reach_ECAL = propagator.getSuccess(); // 0 does not reach ECAL, 1 yes barrel, 2 yes endcaps
127  // ECAL entry point for track
128  GlobalPoint ecal_pos(propagator.particle().x(), propagator.particle().y(), propagator.particle().z());
129 
130  // Track-cluster matching for most energetic clusters
131  sc_clus1_nxtal = -999;
132  sc_clus1_dphi = -999.;
133  sc_clus2_dphi = -999.;
134  sc_clus1_deta = -999.;
135  sc_clus2_deta = -999.;
136  sc_clus1_E = -999.;
137  sc_clus2_E = -999.;
138  sc_clus1_E_ov_p = -999.;
139  sc_clus2_E_ov_p = -999.;
141  *gsf,
142  reach_ECAL,
143  ecal_pos,
144  sc_clus1_nxtal,
145  sc_clus1_dphi,
146  sc_clus2_dphi,
147  sc_clus1_deta,
148  sc_clus2_deta,
149  sc_clus1_E,
150  sc_clus2_E,
151  sc_clus1_E_ov_p,
152  sc_clus2_E_ov_p);
153  sc_clus1_nxtal = (int)sc_clus1_nxtal;
154 
155  } // sc.isNonnull()
156  } // gsf.isNonnull()
157  } // clusters
158 
159  // Out-of-range
160  eid_rho = std::clamp(eid_rho, 0.f, 100.f);
161  eid_sc_eta = std::clamp(eid_sc_eta, -5.f, 5.f);
162  eid_shape_full5x5_r9 = std::clamp(eid_shape_full5x5_r9, 0.f, 2.f);
163  eid_sc_etaWidth = std::clamp(eid_sc_etaWidth, 0.f, 3.14f);
164  eid_sc_phiWidth = std::clamp(eid_sc_phiWidth, 0.f, 3.14f);
165  eid_shape_full5x5_HoverE = std::clamp(eid_shape_full5x5_HoverE, 0.f, 50.f);
166  eid_trk_nhits = std::clamp(eid_trk_nhits, -1.f, 50.f);
167  eid_trk_chi2red = std::clamp(eid_trk_chi2red, -1.f, 50.f);
168  eid_gsf_chi2red = std::clamp(eid_gsf_chi2red, -1.f, 100.f);
169  if (eid_brem_frac < 0.)
170  eid_brem_frac = -1.; //
171  if (eid_brem_frac > 1.)
172  eid_brem_frac = 1.; //
173  eid_gsf_nhits = std::clamp(eid_gsf_nhits, -1.f, 50.f);
174  eid_match_SC_EoverP = std::clamp(eid_match_SC_EoverP, 0.f, 100.f);
175  eid_match_eclu_EoverP = std::clamp(eid_match_eclu_EoverP, -1.f, 1.f);
176  eid_match_SC_dEta = std::clamp(eid_match_SC_dEta, -10.f, 10.f);
177  eid_match_SC_dPhi = std::clamp(eid_match_SC_dPhi, -3.14f, 3.14f);
178  eid_match_seed_dEta = std::clamp(eid_match_seed_dEta, -10.f, 10.f);
179  eid_sc_E = std::clamp(eid_sc_E, 0.f, 1000.f);
180  eid_trk_p = std::clamp(eid_trk_p, -1.f, 1000.f);
181  gsf_mode_p = std::clamp(gsf_mode_p, 0.f, 1000.f);
182  core_shFracHits = std::clamp(core_shFracHits, 0.f, 1.f);
183  gsf_bdtout1 = std::clamp(gsf_bdtout1, -20.f, 20.f);
184  if (gsf_dr < 0.)
185  gsf_dr = 5.; //
186  if (gsf_dr > 5.)
187  gsf_dr = 5.; //
188  if (trk_dr < 0.)
189  trk_dr = 5.; //
190  if (trk_dr > 5.)
191  trk_dr = 5.; //
192  sc_Nclus = std::clamp(sc_Nclus, 0.f, 20.f);
193  sc_clus1_nxtal = std::clamp(sc_clus1_nxtal, 0.f, 100.f);
194  sc_clus1_dphi = std::clamp(sc_clus1_dphi, -3.14f, 3.14f);
195  sc_clus2_dphi = std::clamp(sc_clus2_dphi, -3.14f, 3.14f);
196  sc_clus1_deta = std::clamp(sc_clus1_deta, -5.f, 5.f);
197  sc_clus2_deta = std::clamp(sc_clus2_deta, -5.f, 5.f);
198  sc_clus1_E = std::clamp(sc_clus1_E, 0.f, 1000.f);
199  sc_clus2_E = std::clamp(sc_clus2_E, 0.f, 1000.f);
200  if (sc_clus1_E_ov_p < 0.)
201  sc_clus1_E_ov_p = -1.; //
202  if (sc_clus2_E_ov_p < 0.)
203  sc_clus2_E_ov_p = -1.; //
204 
205  // Set contents of vector
206  std::vector<float> output = {eid_rho,
207  eid_sc_eta,
208  eid_shape_full5x5_r9,
209  eid_sc_etaWidth,
210  eid_sc_phiWidth,
211  eid_shape_full5x5_HoverE,
212  eid_trk_nhits,
213  eid_trk_chi2red,
214  eid_gsf_chi2red,
215  eid_brem_frac,
216  eid_gsf_nhits,
217  eid_match_SC_EoverP,
218  eid_match_eclu_EoverP,
219  eid_match_SC_dEta,
220  eid_match_SC_dPhi,
221  eid_match_seed_dEta,
222  eid_sc_E,
223  eid_trk_p,
224  gsf_mode_p,
225  core_shFracHits,
226  gsf_bdtout1,
227  gsf_dr,
228  trk_dr,
229  sc_Nclus,
230  sc_clus1_nxtal,
231  sc_clus1_dphi,
232  sc_clus2_dphi,
233  sc_clus1_deta,
234  sc_clus2_deta,
235  sc_clus1_E,
236  sc_clus2_E,
237  sc_clus1_E_ov_p,
238  sc_clus2_E_ov_p};
239  return output;
240  }
241 
243  // feature list for original models (2019Aug07 and earlier)
244  std::vector<float> features_V0(reco::GsfElectron const& ele, float rho, float unbiased) {
245  float eid_rho = -999.;
246  float eid_sc_eta = -999.;
247  float eid_shape_full5x5_r9 = -999.;
248  float eid_sc_etaWidth = -999.;
249  float eid_sc_phiWidth = -999.;
250  float eid_shape_full5x5_HoverE = -999.;
251  float eid_trk_nhits = -999.;
252  float eid_trk_chi2red = -999.;
253  float eid_gsf_chi2red = -999.;
254  float eid_brem_frac = -999.;
255  float eid_gsf_nhits = -999.;
256  float eid_match_SC_EoverP = -999.;
257  float eid_match_eclu_EoverP = -999.;
258  float eid_match_SC_dEta = -999.;
259  float eid_match_SC_dPhi = -999.;
260  float eid_match_seed_dEta = -999.;
261  float eid_sc_E = -999.;
262  float eid_trk_p = -999.;
263  float gsf_mode_p = -999.;
264  float core_shFracHits = -999.;
265  float gsf_bdtout1 = -999.;
266  float gsf_dr = -999.;
267  float trk_dr = -999.;
268  float sc_Nclus = -999.;
269  float sc_clus1_nxtal = -999.;
270  float sc_clus1_dphi = -999.;
271  float sc_clus2_dphi = -999.;
272  float sc_clus1_deta = -999.;
273  float sc_clus2_deta = -999.;
274  float sc_clus1_E = -999.;
275  float sc_clus2_E = -999.;
276  float sc_clus1_E_ov_p = -999.;
277  float sc_clus2_E_ov_p = -999.;
278 
279  // KF tracks
280  if (ele.core().isNonnull()) {
281  const auto& trk = ele.closestCtfTrackRef(); // reco::TrackRef
282  if (trk.isNonnull()) {
283  eid_trk_p = (float)trk->p();
284  eid_trk_nhits = (float)trk->found();
285  eid_trk_chi2red = (float)trk->normalizedChi2();
286  TVector3 trkTV3(0, 0, 0);
287  trkTV3.SetPtEtaPhi(trk->pt(), trk->eta(), trk->phi());
288  TVector3 eleTV3(0, 0, 0);
289  eleTV3.SetPtEtaPhi(ele.pt(), ele.eta(), ele.phi());
290  trk_dr = eleTV3.DeltaR(trkTV3);
291  }
292  }
293 
294  // GSF tracks
295  if (ele.core().isNonnull()) {
296  const auto& gsf = ele.core()->gsfTrack(); // reco::GsfTrackRef
297  if (gsf.isNonnull()) {
298  gsf_mode_p = gsf->pMode();
299  eid_gsf_nhits = (float)gsf->found();
300  eid_gsf_chi2red = gsf->normalizedChi2();
301  TVector3 gsfTV3(0, 0, 0);
302  gsfTV3.SetPtEtaPhi(gsf->ptMode(), gsf->etaMode(), gsf->phiMode());
303  TVector3 eleTV3(0, 0, 0);
304  eleTV3.SetPtEtaPhi(ele.pt(), ele.eta(), ele.phi());
305  gsf_dr = eleTV3.DeltaR(gsfTV3);
306  }
307  }
308 
309  // Super clusters
310  if (ele.core().isNonnull()) {
311  const auto& sc = ele.core()->superCluster(); // reco::SuperClusterRef
312  if (sc.isNonnull()) {
313  eid_sc_E = sc->energy();
314  eid_sc_eta = sc->eta();
315  eid_sc_etaWidth = sc->etaWidth();
316  eid_sc_phiWidth = sc->phiWidth();
317  sc_Nclus = (float)sc->clustersSize();
318  }
319  }
320 
321  // Track-cluster matching
322  eid_match_seed_dEta = ele.deltaEtaSeedClusterTrackAtCalo();
323  eid_match_eclu_EoverP = (1. / ele.ecalEnergy()) - (1. / ele.p());
324  eid_match_SC_EoverP = ele.eSuperClusterOverP();
325  eid_match_SC_dEta = ele.deltaEtaSuperClusterTrackAtVtx();
326  eid_match_SC_dPhi = ele.deltaPhiSuperClusterTrackAtVtx();
327 
328  // Shower shape vars
329  eid_shape_full5x5_HoverE = ele.full5x5_hcalOverEcal();
330  eid_shape_full5x5_r9 = ele.full5x5_r9();
331 
332  // Misc
333  eid_rho = rho;
334 
335  eid_brem_frac = ele.fbrem();
336  core_shFracHits = (float)ele.shFracInnerHits();
337 
338  // Unbiased BDT from ElectronSeed
339  gsf_bdtout1 = unbiased;
340 
341  // Clusters
342  if (ele.core().isNonnull()) {
343  const auto& gsf = ele.core()->gsfTrack(); // reco::GsfTrackRef
344  if (gsf.isNonnull()) {
345  const auto& sc = ele.core()->superCluster(); // reco::SuperClusterRef
346  if (sc.isNonnull()) {
347  // Propagate electron track to ECAL surface
348  double mass2 = 0.000511 * 0.000511;
349  float p2 = pow(gsf->p(), 2);
350  float energy = sqrt(mass2 + p2);
351  math::XYZTLorentzVector mom = math::XYZTLorentzVector(gsf->px(), gsf->py(), gsf->pz(), energy);
352  math::XYZTLorentzVector pos = math::XYZTLorentzVector(gsf->vx(), gsf->vy(), gsf->vz(), 0.);
353  float field_z = 3.8;
354  BaseParticlePropagator mypart(RawParticle(mom, pos, gsf->charge()), 0, 0, field_z);
355  mypart.propagateToEcalEntrance(true); // true only first half loop , false more than one loop
356  bool reach_ECAL = mypart.getSuccess(); // 0 does not reach ECAL, 1 yes barrel, 2 yes endcaps
357 
358  // ECAL entry point for track
359  GlobalPoint ecal_pos(
360  mypart.particle().vertex().x(), mypart.particle().vertex().y(), mypart.particle().vertex().z());
361 
362  // Track-cluster matching for most energetic clusters
363  sc_clus1_nxtal = -999.;
364  sc_clus1_dphi = -999.;
365  sc_clus2_dphi = -999.;
366  sc_clus1_deta = -999.;
367  sc_clus2_deta = -999.;
368  sc_clus1_E = -999.;
369  sc_clus2_E = -999.;
370  sc_clus1_E_ov_p = -999.;
371  sc_clus2_E_ov_p = -999.;
373  *gsf,
374  reach_ECAL,
375  ecal_pos,
376  sc_clus1_nxtal,
377  sc_clus1_dphi,
378  sc_clus2_dphi,
379  sc_clus1_deta,
380  sc_clus2_deta,
381  sc_clus1_E,
382  sc_clus2_E,
383  sc_clus1_E_ov_p,
384  sc_clus2_E_ov_p);
385 
386  } // sc.isNonnull()
387  } // gsf.isNonnull()
388  } // clusters
389 
390  // Out-of-range
391  eid_sc_eta = std::clamp(eid_sc_eta, -5.f, 5.f);
392  eid_shape_full5x5_r9 = std::clamp(eid_shape_full5x5_r9, 0.f, 2.f);
393  eid_sc_etaWidth = std::clamp(eid_sc_etaWidth, 0.f, 3.14f);
394  eid_sc_phiWidth = std::clamp(eid_sc_phiWidth, 0.f, 3.14f);
395  eid_shape_full5x5_HoverE = std::clamp(eid_shape_full5x5_HoverE, 0.f, 50.f);
396  eid_trk_nhits = std::clamp(eid_trk_nhits, -1.f, 50.f);
397  eid_trk_chi2red = std::clamp(eid_trk_chi2red, -1.f, 50.f);
398  eid_gsf_chi2red = std::clamp(eid_gsf_chi2red, -1.f, 100.f);
399  if (eid_brem_frac < 0.)
400  eid_brem_frac = -1.; //
401  if (eid_brem_frac > 1.)
402  eid_brem_frac = 1.; //
403  eid_gsf_nhits = std::clamp(eid_gsf_nhits, -1.f, 50.f);
404  eid_match_SC_EoverP = std::clamp(eid_match_SC_EoverP, 0.f, 100.f);
405  eid_match_eclu_EoverP = std::clamp(eid_match_eclu_EoverP, -1.f, 1.f);
406  eid_match_SC_dEta = std::clamp(eid_match_SC_dEta, -10.f, 10.f);
407  eid_match_SC_dPhi = std::clamp(eid_match_SC_dPhi, -3.14f, 3.14f);
408  eid_match_seed_dEta = std::clamp(eid_match_seed_dEta, -10.f, 10.f);
409  eid_sc_E = std::clamp(eid_sc_E, 0.f, 1000.f);
410  eid_trk_p = std::clamp(eid_trk_p, -1.f, 1000.f);
411  gsf_mode_p = std::clamp(gsf_mode_p, 0.f, 1000.f);
412  core_shFracHits = std::clamp(core_shFracHits, 0.f, 1.f);
413  gsf_bdtout1 = std::clamp(gsf_bdtout1, -20.f, 20.f);
414  if (gsf_dr < 0.)
415  gsf_dr = 5.; //
416  if (gsf_dr > 5.)
417  gsf_dr = 5.; //
418  if (trk_dr < 0.)
419  trk_dr = 5.; //
420  if (trk_dr > 5.)
421  trk_dr = 5.; //
422  sc_Nclus = std::clamp(sc_Nclus, 0.f, 20.f);
423  sc_clus1_nxtal = std::clamp(sc_clus1_nxtal, 0.f, 100.f);
424  if (sc_clus1_dphi < -3.14)
425  sc_clus1_dphi = -5.; //
426  if (sc_clus1_dphi > 3.14)
427  sc_clus1_dphi = 5.; //
428  if (sc_clus2_dphi < -3.14)
429  sc_clus2_dphi = -5.; //
430  if (sc_clus2_dphi > 3.14)
431  sc_clus2_dphi = 5.; //
432  sc_clus1_deta = std::clamp(sc_clus1_deta, -5.f, 5.f);
433  sc_clus2_deta = std::clamp(sc_clus2_deta, -5.f, 5.f);
434  sc_clus1_E = std::clamp(sc_clus1_E, 0.f, 1000.f);
435  sc_clus2_E = std::clamp(sc_clus2_E, 0.f, 1000.f);
436  if (sc_clus1_E_ov_p < 0.)
437  sc_clus1_E_ov_p = -1.; //
438  if (sc_clus2_E_ov_p < 0.)
439  sc_clus2_E_ov_p = -1.; //
440 
441  // Set contents of vector
442  std::vector<float> output = {eid_rho,
443  eid_sc_eta,
444  eid_shape_full5x5_r9,
445  eid_sc_etaWidth,
446  eid_sc_phiWidth,
447  eid_shape_full5x5_HoverE,
448  eid_trk_nhits,
449  eid_trk_chi2red,
450  eid_gsf_chi2red,
451  eid_brem_frac,
452  eid_gsf_nhits,
453  eid_match_SC_EoverP,
454  eid_match_eclu_EoverP,
455  eid_match_SC_dEta,
456  eid_match_SC_dPhi,
457  eid_match_seed_dEta,
458  eid_sc_E,
459  eid_trk_p,
460  gsf_mode_p,
461  core_shFracHits,
462  gsf_bdtout1,
463  gsf_dr,
464  trk_dr,
465  sc_Nclus,
466  sc_clus1_nxtal,
467  sc_clus1_dphi,
468  sc_clus2_dphi,
469  sc_clus1_deta,
470  sc_clus2_deta,
471  sc_clus1_E,
472  sc_clus2_E,
473  sc_clus1_E_ov_p,
474  sc_clus2_E_ov_p};
475  return output;
476  }
477 
479  // Find most energetic clusters
481  reco::SuperCluster const& sc, int& clusNum, float& maxEne1, float& maxEne2, int& i1, int& i2) {
482  if (sc.clustersSize() > 0 && sc.clustersBegin() != sc.clustersEnd()) {
483  for (auto const& cluster : sc.clusters()) {
484  if (cluster->energy() > maxEne1) {
485  maxEne1 = cluster->energy();
486  i1 = clusNum;
487  }
488  clusNum++;
489  }
490  if (sc.clustersSize() > 1) {
491  clusNum = 0;
492  for (auto const& cluster : sc.clusters()) {
493  if (clusNum != i1) {
494  if (cluster->energy() > maxEne2) {
495  maxEne2 = cluster->energy();
496  i2 = clusNum;
497  }
498  }
499  clusNum++;
500  }
501  }
502  } // loop over clusters
503  }
504 
506  // Track-cluster matching for most energetic clusters
508  reco::GsfTrack const& gsf,
509  bool const& reach_ECAL,
510  GlobalPoint const& ecal_pos,
511  float& sc_clus1_nxtal,
512  float& sc_clus1_dphi,
513  float& sc_clus2_dphi,
514  float& sc_clus1_deta,
515  float& sc_clus2_deta,
516  float& sc_clus1_E,
517  float& sc_clus2_E,
518  float& sc_clus1_E_ov_p,
519  float& sc_clus2_E_ov_p) {
520  // Iterate through ECAL clusters and sort in energy
521  int clusNum = 0;
522  float maxEne1 = -1;
523  float maxEne2 = -1;
524  int i1 = -1;
525  int i2 = -1;
526  findEnergeticClusters(sc, clusNum, maxEne1, maxEne2, i1, i2);
527 
528  // track-clusters match
529  clusNum = 0;
530  if (sc.clustersSize() > 0 && sc.clustersBegin() != sc.clustersEnd()) {
531  for (auto const& cluster : sc.clusters()) {
532  float deta = ecal_pos.eta() - cluster->eta();
533  float dphi = reco::deltaPhi(ecal_pos.phi(), cluster->phi());
534  if (clusNum == i1) {
535  sc_clus1_E = cluster->energy();
536  if (gsf.pMode() > 0)
537  sc_clus1_E_ov_p = cluster->energy() / gsf.pMode();
538  sc_clus1_nxtal = (float)cluster->size();
539  if (reach_ECAL > 0) {
540  sc_clus1_deta = deta;
541  sc_clus1_dphi = dphi;
542  }
543  } else if (clusNum == i2) {
544  sc_clus2_E = cluster->energy();
545  if (gsf.pMode() > 0)
546  sc_clus2_E_ov_p = cluster->energy() / gsf.pMode();
547  if (reach_ECAL > 0) {
548  sc_clus2_deta = deta;
549  sc_clus2_dphi = dphi;
550  }
551  }
552  clusNum++;
553  }
554  }
555  }
556 
557 } // namespace lowptgsfeleid
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:186
virtual TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:205
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
void findEnergeticClusters(reco::SuperCluster const &, int &, float &, float &, int &, int &)
double eta() const final
momentum pseudorapidity
double pMode() const
momentum vector magnitude from mode
Definition: GsfTrack.h:47
float eSuperClusterOverP() const
Definition: GsfElectron.h:249
double y() const
y of vertex
Definition: RawParticle.h:302
double z() const
z of vertex
Definition: RawParticle.h:303
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
int getSuccess() const
Has propagation been performed and was barrel or endcap reached ?
double pt() const final
transverse momentum
float fbrem() const
Definition: GsfElectron.h:772
RawParticle const & particle() const
The particle being propagated.
const CaloClusterPtrVector & clusters() const
const access to the cluster list itself
Definition: SuperCluster.h:69
void trackClusterMatching(reco::SuperCluster const &, reco::GsfTrack const &, bool const &, GlobalPoint const &, float &, float &, float &, float &, float &, float &, float &, float &, float &)
double x() const
x of vertex
Definition: RawParticle.h:301
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:253
T sqrt(T t)
Definition: SSEVec.h:18
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:256
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
double f[11][100]
std::vector< float > features_V0(reco::GsfElectron const &ele, float rho, float unbiased)
double p2[4]
Definition: TauolaWrapper.h:90
bool propagateToEcalEntrance(bool first=true)
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:339
float shFracInnerHits() const
Definition: GsfElectron.h:204
double p() const final
magnitude of momentum vector
size_t clustersSize() const
number of BasicCluster constituents
Definition: SuperCluster.h:87
float full5x5_hcalOverEcal() const
Definition: GsfElectron.h:470
float ecalEnergy() const
Definition: GsfElectron.h:860
virtual GsfElectronCoreRef core() const
Definition: GsfElectron.cc:8
float full5x5_r9() const
Definition: GsfElectron.h:467
T eta() const
Definition: PV3DBase.h:76
float deltaEtaSeedClusterTrackAtCalo() const
Definition: GsfElectron.h:254
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:185
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:75
double phi() const final
momentum azimuthal angle
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::vector< float > features_V1(reco::GsfElectron const &ele, float rho, float unbiased, float field_z)
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:78
float eta() const
Momentum pseudorapidity. Note this is taken from the simtrack before the calorimeter.
Definition: SimCluster.h:148