CMS 3D CMS Logo

DeepTauId.cc
Go to the documentation of this file.
1 /*
2  * \class DeepTauId
3  *
4  * Tau identification using Deep NN.
5  *
6  * \author Konstantin Androsov, INFN Pisa
7  */
8 
11 
12 namespace deep_tau {
13  constexpr int NumberOfOutputs = 4;
14 }
15 
16 namespace {
17 
18  struct dnn_inputs_2017v1 {
19  enum vars {
20  pt = 0,
21  eta,
22  mass,
23  decayMode,
24  chargedIsoPtSum,
25  neutralIsoPtSum,
26  neutralIsoPtSumWeight,
27  photonPtSumOutsideSignalCone,
28  puCorrPtSum,
29  dxy,
30  dxy_sig,
31  dz,
32  ip3d,
33  ip3d_sig,
34  hasSecondaryVertex,
35  flightLength_r,
36  flightLength_dEta,
37  flightLength_dPhi,
38  flightLength_sig,
39  leadChargedHadrCand_pt,
40  leadChargedHadrCand_dEta,
41  leadChargedHadrCand_dPhi,
42  leadChargedHadrCand_mass,
47  leadingTrackNormChi2,
48  e_ratio,
49  gj_angle_diff,
50  n_photons,
51  emFraction,
52  has_gsf_track,
53  inside_ecal_crack,
54  gsf_ele_matched,
55  gsf_ele_pt,
56  gsf_ele_dEta,
57  gsf_ele_dPhi,
58  gsf_ele_mass,
59  gsf_ele_Ee,
60  gsf_ele_Egamma,
61  gsf_ele_Pin,
62  gsf_ele_Pout,
63  gsf_ele_EtotOverPin,
64  gsf_ele_Eecal,
65  gsf_ele_dEta_SeedClusterTrackAtCalo,
66  gsf_ele_dPhi_SeedClusterTrackAtCalo,
67  gsf_ele_mvaIn_sigmaEtaEta,
68  gsf_ele_mvaIn_hadEnergy,
69  gsf_ele_mvaIn_deltaEta,
70  gsf_ele_Chi2NormGSF,
71  gsf_ele_GSFNumHits,
72  gsf_ele_GSFTrackResol,
73  gsf_ele_GSFTracklnPt,
74  gsf_ele_Chi2NormKF,
75  gsf_ele_KFNumHits,
76  leadChargedCand_etaAtEcalEntrance,
77  leadChargedCand_pt,
78  leadChargedHadrCand_HoP,
79  leadChargedHadrCand_EoP,
80  tau_visMass_innerSigCone,
81  n_matched_muons,
82  muon_pt,
83  muon_dEta,
84  muon_dPhi,
85  muon_n_matches_DT_1,
86  muon_n_matches_DT_2,
87  muon_n_matches_DT_3,
88  muon_n_matches_DT_4,
89  muon_n_matches_CSC_1,
90  muon_n_matches_CSC_2,
91  muon_n_matches_CSC_3,
92  muon_n_matches_CSC_4,
93  muon_n_hits_DT_2,
94  muon_n_hits_DT_3,
95  muon_n_hits_DT_4,
96  muon_n_hits_CSC_2,
97  muon_n_hits_CSC_3,
98  muon_n_hits_CSC_4,
99  muon_n_hits_RPC_2,
100  muon_n_hits_RPC_3,
101  muon_n_hits_RPC_4,
102  muon_n_stations_with_matches_03,
103  muon_n_stations_with_hits_23,
104  signalChargedHadrCands_sum_innerSigCone_pt,
105  signalChargedHadrCands_sum_innerSigCone_dEta,
106  signalChargedHadrCands_sum_innerSigCone_dPhi,
107  signalChargedHadrCands_sum_innerSigCone_mass,
108  signalChargedHadrCands_sum_outerSigCone_pt,
109  signalChargedHadrCands_sum_outerSigCone_dEta,
110  signalChargedHadrCands_sum_outerSigCone_dPhi,
111  signalChargedHadrCands_sum_outerSigCone_mass,
112  signalChargedHadrCands_nTotal_innerSigCone,
113  signalChargedHadrCands_nTotal_outerSigCone,
114  signalNeutrHadrCands_sum_innerSigCone_pt,
115  signalNeutrHadrCands_sum_innerSigCone_dEta,
116  signalNeutrHadrCands_sum_innerSigCone_dPhi,
117  signalNeutrHadrCands_sum_innerSigCone_mass,
118  signalNeutrHadrCands_sum_outerSigCone_pt,
119  signalNeutrHadrCands_sum_outerSigCone_dEta,
120  signalNeutrHadrCands_sum_outerSigCone_dPhi,
121  signalNeutrHadrCands_sum_outerSigCone_mass,
122  signalNeutrHadrCands_nTotal_innerSigCone,
123  signalNeutrHadrCands_nTotal_outerSigCone,
124  signalGammaCands_sum_innerSigCone_pt,
125  signalGammaCands_sum_innerSigCone_dEta,
126  signalGammaCands_sum_innerSigCone_dPhi,
127  signalGammaCands_sum_innerSigCone_mass,
128  signalGammaCands_sum_outerSigCone_pt,
129  signalGammaCands_sum_outerSigCone_dEta,
130  signalGammaCands_sum_outerSigCone_dPhi,
131  signalGammaCands_sum_outerSigCone_mass,
132  signalGammaCands_nTotal_innerSigCone,
133  signalGammaCands_nTotal_outerSigCone,
134  isolationChargedHadrCands_sum_pt,
135  isolationChargedHadrCands_sum_dEta,
136  isolationChargedHadrCands_sum_dPhi,
137  isolationChargedHadrCands_sum_mass,
138  isolationChargedHadrCands_nTotal,
139  isolationNeutrHadrCands_sum_pt,
140  isolationNeutrHadrCands_sum_dEta,
141  isolationNeutrHadrCands_sum_dPhi,
142  isolationNeutrHadrCands_sum_mass,
143  isolationNeutrHadrCands_nTotal,
144  isolationGammaCands_sum_pt,
145  isolationGammaCands_sum_dEta,
146  isolationGammaCands_sum_dPhi,
147  isolationGammaCands_sum_mass,
148  isolationGammaCands_nTotal,
149  NumberOfInputs
150  };
151  };
152 
153  namespace dnn_inputs_2017_v2 {
154  constexpr int number_of_inner_cell = 11;
155  constexpr int number_of_outer_cell = 21;
156  constexpr int number_of_conv_features = 64;
157  namespace TauBlockInputs {
158  enum vars {
159  rho = 0,
160  tau_pt,
161  tau_eta,
162  tau_phi,
163  tau_mass,
164  tau_E_over_pt,
165  tau_charge,
166  tau_n_charged_prongs,
167  tau_n_neutral_prongs,
168  chargedIsoPtSum,
169  chargedIsoPtSumdR03_over_dR05,
170  footprintCorrection,
171  neutralIsoPtSum,
172  neutralIsoPtSumWeight_over_neutralIsoPtSum,
173  neutralIsoPtSumWeightdR03_over_neutralIsoPtSum,
174  neutralIsoPtSumdR03_over_dR05,
175  photonPtSumOutsideSignalCone,
176  puCorrPtSum,
177  tau_dxy_pca_x,
178  tau_dxy_pca_y,
179  tau_dxy_pca_z,
180  tau_dxy_valid,
181  tau_dxy,
182  tau_dxy_sig,
183  tau_ip3d_valid,
184  tau_ip3d,
185  tau_ip3d_sig,
186  tau_dz,
187  tau_dz_sig_valid,
188  tau_dz_sig,
189  tau_flightLength_x,
190  tau_flightLength_y,
191  tau_flightLength_z,
192  tau_flightLength_sig,
193  tau_pt_weighted_deta_strip,
194  tau_pt_weighted_dphi_strip,
195  tau_pt_weighted_dr_signal,
196  tau_pt_weighted_dr_iso,
197  tau_leadingTrackNormChi2,
198  tau_e_ratio_valid,
199  tau_e_ratio,
200  tau_gj_angle_diff_valid,
201  tau_gj_angle_diff,
202  tau_n_photons,
203  tau_emFraction,
204  tau_inside_ecal_crack,
205  leadChargedCand_etaAtEcalEntrance_minus_tau_eta,
206  NumberOfInputs
207  };
208  }
209 
210  namespace EgammaBlockInputs {
211  enum vars {
212  rho = 0,
213  tau_pt,
214  tau_eta,
215  tau_inside_ecal_crack,
216  pfCand_ele_valid,
217  pfCand_ele_rel_pt,
218  pfCand_ele_deta,
219  pfCand_ele_dphi,
220  pfCand_ele_pvAssociationQuality,
221  pfCand_ele_puppiWeight,
222  pfCand_ele_charge,
223  pfCand_ele_lostInnerHits,
224  pfCand_ele_numberOfPixelHits,
225  pfCand_ele_vertex_dx,
226  pfCand_ele_vertex_dy,
227  pfCand_ele_vertex_dz,
228  pfCand_ele_vertex_dx_tauFL,
229  pfCand_ele_vertex_dy_tauFL,
230  pfCand_ele_vertex_dz_tauFL,
231  pfCand_ele_hasTrackDetails,
232  pfCand_ele_dxy,
233  pfCand_ele_dxy_sig,
234  pfCand_ele_dz,
235  pfCand_ele_dz_sig,
236  pfCand_ele_track_chi2_ndof,
237  pfCand_ele_track_ndof,
238  ele_valid,
239  ele_rel_pt,
240  ele_deta,
241  ele_dphi,
242  ele_cc_valid,
243  ele_cc_ele_rel_energy,
244  ele_cc_gamma_rel_energy,
245  ele_cc_n_gamma,
246  ele_rel_trackMomentumAtVtx,
247  ele_rel_trackMomentumAtCalo,
248  ele_rel_trackMomentumOut,
249  ele_rel_trackMomentumAtEleClus,
250  ele_rel_trackMomentumAtVtxWithConstraint,
251  ele_rel_ecalEnergy,
252  ele_ecalEnergy_sig,
253  ele_eSuperClusterOverP,
254  ele_eSeedClusterOverP,
255  ele_eSeedClusterOverPout,
256  ele_eEleClusterOverPout,
257  ele_deltaEtaSuperClusterTrackAtVtx,
258  ele_deltaEtaSeedClusterTrackAtCalo,
259  ele_deltaEtaEleClusterTrackAtCalo,
260  ele_deltaPhiEleClusterTrackAtCalo,
261  ele_deltaPhiSuperClusterTrackAtVtx,
262  ele_deltaPhiSeedClusterTrackAtCalo,
263  ele_mvaInput_earlyBrem,
264  ele_mvaInput_lateBrem,
265  ele_mvaInput_sigmaEtaEta,
266  ele_mvaInput_hadEnergy,
267  ele_mvaInput_deltaEta,
268  ele_gsfTrack_normalizedChi2,
269  ele_gsfTrack_numberOfValidHits,
270  ele_rel_gsfTrack_pt,
271  ele_gsfTrack_pt_sig,
272  ele_has_closestCtfTrack,
273  ele_closestCtfTrack_normalizedChi2,
274  ele_closestCtfTrack_numberOfValidHits,
275  pfCand_gamma_valid,
276  pfCand_gamma_rel_pt,
277  pfCand_gamma_deta,
278  pfCand_gamma_dphi,
279  pfCand_gamma_pvAssociationQuality,
280  pfCand_gamma_fromPV,
281  pfCand_gamma_puppiWeight,
282  pfCand_gamma_puppiWeightNoLep,
283  pfCand_gamma_lostInnerHits,
284  pfCand_gamma_numberOfPixelHits,
285  pfCand_gamma_vertex_dx,
286  pfCand_gamma_vertex_dy,
287  pfCand_gamma_vertex_dz,
288  pfCand_gamma_vertex_dx_tauFL,
289  pfCand_gamma_vertex_dy_tauFL,
290  pfCand_gamma_vertex_dz_tauFL,
291  pfCand_gamma_hasTrackDetails,
292  pfCand_gamma_dxy,
293  pfCand_gamma_dxy_sig,
294  pfCand_gamma_dz,
295  pfCand_gamma_dz_sig,
296  pfCand_gamma_track_chi2_ndof,
297  pfCand_gamma_track_ndof,
298  NumberOfInputs
299  };
300  }
301 
302  namespace MuonBlockInputs {
303  enum vars {
304  rho = 0,
305  tau_pt,
306  tau_eta,
307  tau_inside_ecal_crack,
308  pfCand_muon_valid,
309  pfCand_muon_rel_pt,
310  pfCand_muon_deta,
311  pfCand_muon_dphi,
312  pfCand_muon_pvAssociationQuality,
313  pfCand_muon_fromPV,
314  pfCand_muon_puppiWeight,
315  pfCand_muon_charge,
316  pfCand_muon_lostInnerHits,
317  pfCand_muon_numberOfPixelHits,
318  pfCand_muon_vertex_dx,
319  pfCand_muon_vertex_dy,
320  pfCand_muon_vertex_dz,
321  pfCand_muon_vertex_dx_tauFL,
322  pfCand_muon_vertex_dy_tauFL,
323  pfCand_muon_vertex_dz_tauFL,
324  pfCand_muon_hasTrackDetails,
325  pfCand_muon_dxy,
326  pfCand_muon_dxy_sig,
327  pfCand_muon_dz,
328  pfCand_muon_dz_sig,
329  pfCand_muon_track_chi2_ndof,
330  pfCand_muon_track_ndof,
331  muon_valid,
332  muon_rel_pt,
333  muon_deta,
334  muon_dphi,
335  muon_dxy,
336  muon_dxy_sig,
337  muon_normalizedChi2_valid,
338  muon_normalizedChi2,
339  muon_numberOfValidHits,
340  muon_segmentCompatibility,
341  muon_caloCompatibility,
342  muon_pfEcalEnergy_valid,
343  muon_rel_pfEcalEnergy,
344  muon_n_matches_DT_1,
345  muon_n_matches_DT_2,
346  muon_n_matches_DT_3,
347  muon_n_matches_DT_4,
348  muon_n_matches_CSC_1,
349  muon_n_matches_CSC_2,
350  muon_n_matches_CSC_3,
351  muon_n_matches_CSC_4,
352  muon_n_matches_RPC_1,
353  muon_n_matches_RPC_2,
354  muon_n_matches_RPC_3,
355  muon_n_matches_RPC_4,
356  muon_n_hits_DT_1,
357  muon_n_hits_DT_2,
358  muon_n_hits_DT_3,
359  muon_n_hits_DT_4,
360  muon_n_hits_CSC_1,
361  muon_n_hits_CSC_2,
362  muon_n_hits_CSC_3,
363  muon_n_hits_CSC_4,
364  muon_n_hits_RPC_1,
365  muon_n_hits_RPC_2,
366  muon_n_hits_RPC_3,
367  muon_n_hits_RPC_4,
368  NumberOfInputs
369  };
370  }
371 
372  namespace HadronBlockInputs {
373  enum vars {
374  rho = 0,
375  tau_pt,
376  tau_eta,
377  tau_inside_ecal_crack,
378  pfCand_chHad_valid,
379  pfCand_chHad_rel_pt,
380  pfCand_chHad_deta,
381  pfCand_chHad_dphi,
382  pfCand_chHad_leadChargedHadrCand,
383  pfCand_chHad_pvAssociationQuality,
384  pfCand_chHad_fromPV,
385  pfCand_chHad_puppiWeight,
386  pfCand_chHad_puppiWeightNoLep,
387  pfCand_chHad_charge,
388  pfCand_chHad_lostInnerHits,
389  pfCand_chHad_numberOfPixelHits,
390  pfCand_chHad_vertex_dx,
391  pfCand_chHad_vertex_dy,
392  pfCand_chHad_vertex_dz,
393  pfCand_chHad_vertex_dx_tauFL,
394  pfCand_chHad_vertex_dy_tauFL,
395  pfCand_chHad_vertex_dz_tauFL,
396  pfCand_chHad_hasTrackDetails,
397  pfCand_chHad_dxy,
398  pfCand_chHad_dxy_sig,
399  pfCand_chHad_dz,
400  pfCand_chHad_dz_sig,
401  pfCand_chHad_track_chi2_ndof,
402  pfCand_chHad_track_ndof,
403  pfCand_chHad_hcalFraction,
404  pfCand_chHad_rawCaloFraction,
405  pfCand_nHad_valid,
406  pfCand_nHad_rel_pt,
407  pfCand_nHad_deta,
408  pfCand_nHad_dphi,
409  pfCand_nHad_puppiWeight,
410  pfCand_nHad_puppiWeightNoLep,
411  pfCand_nHad_hcalFraction,
412  NumberOfInputs
413  };
414  }
415  } // namespace dnn_inputs_2017_v2
416 
417  template <typename LVector1, typename LVector2>
418  float dEta(const LVector1& p4, const LVector2& tau_p4) {
419  return static_cast<float>(p4.eta() - tau_p4.eta());
420  }
421 
422  template <typename LVector1, typename LVector2>
423  float dPhi(const LVector1& p4_1, const LVector2& p4_2) {
424  return static_cast<float>(reco::deltaPhi(p4_2.phi(), p4_1.phi()));
425  }
426 
427  struct MuonHitMatchV1 {
428  static constexpr int n_muon_stations = 4;
429 
430  std::map<int, std::vector<UInt_t>> n_matches, n_hits;
431  unsigned n_muons{0};
432  const pat::Muon* best_matched_muon{nullptr};
433  double deltaR2_best_match{-1};
434 
435  MuonHitMatchV1() {
436  n_matches[MuonSubdetId::DT].assign(n_muon_stations, 0);
437  n_matches[MuonSubdetId::CSC].assign(n_muon_stations, 0);
438  n_matches[MuonSubdetId::RPC].assign(n_muon_stations, 0);
439  n_hits[MuonSubdetId::DT].assign(n_muon_stations, 0);
440  n_hits[MuonSubdetId::CSC].assign(n_muon_stations, 0);
441  n_hits[MuonSubdetId::RPC].assign(n_muon_stations, 0);
442  }
443 
444  void addMatchedMuon(const pat::Muon& muon, const pat::Tau& tau) {
445  static constexpr int n_stations = 4;
446 
447  ++n_muons;
448  const double dR2 = reco::deltaR2(tau.p4(), muon.p4());
449  if (!best_matched_muon || dR2 < deltaR2_best_match) {
450  best_matched_muon = &muon;
451  deltaR2_best_match = dR2;
452  }
453 
454  for (const auto& segment : muon.matches()) {
455  if (segment.segmentMatches.empty())
456  continue;
457  if (n_matches.count(segment.detector()))
458  ++n_matches.at(segment.detector()).at(segment.station() - 1);
459  }
460 
461  if (muon.outerTrack().isNonnull()) {
462  const auto& hit_pattern = muon.outerTrack()->hitPattern();
463  for (int hit_index = 0; hit_index < hit_pattern.numberOfAllHits(reco::HitPattern::TRACK_HITS); ++hit_index) {
464  auto hit_id = hit_pattern.getHitPattern(reco::HitPattern::TRACK_HITS, hit_index);
465  if (hit_id == 0)
466  break;
467  if (hit_pattern.muonHitFilter(hit_id) && (hit_pattern.getHitType(hit_id) == TrackingRecHit::valid ||
468  hit_pattern.getHitType(hit_id == TrackingRecHit::bad))) {
469  const int station = hit_pattern.getMuonStation(hit_id) - 1;
470  if (station > 0 && station < n_stations) {
471  std::vector<UInt_t>* muon_n_hits = nullptr;
472  if (hit_pattern.muonDTHitFilter(hit_id))
473  muon_n_hits = &n_hits.at(MuonSubdetId::DT);
474  else if (hit_pattern.muonCSCHitFilter(hit_id))
475  muon_n_hits = &n_hits.at(MuonSubdetId::CSC);
476  else if (hit_pattern.muonRPCHitFilter(hit_id))
477  muon_n_hits = &n_hits.at(MuonSubdetId::RPC);
478 
479  if (muon_n_hits)
480  ++muon_n_hits->at(station);
481  }
482  }
483  }
484  }
485  }
486 
487  static std::vector<const pat::Muon*> findMatchedMuons(const pat::Tau& tau,
488  const pat::MuonCollection& muons,
489  double deltaR,
490  double minPt) {
491  const reco::Muon* hadr_cand_muon = nullptr;
492  if (tau.leadPFChargedHadrCand().isNonnull() && tau.leadPFChargedHadrCand()->muonRef().isNonnull())
493  hadr_cand_muon = tau.leadPFChargedHadrCand()->muonRef().get();
494  std::vector<const pat::Muon*> matched_muons;
495  const double dR2 = deltaR * deltaR;
496  for (const pat::Muon& muon : muons) {
497  const reco::Muon* reco_muon = &muon;
498  if (muon.pt() <= minPt)
499  continue;
500  if (reco_muon == hadr_cand_muon)
501  continue;
502  if (reco::deltaR2(tau.p4(), muon.p4()) >= dR2)
503  continue;
504  matched_muons.push_back(&muon);
505  }
506  return matched_muons;
507  }
508 
509  template <typename dnn, typename TensorElemGet>
510  void fillTensor(const TensorElemGet& get, const pat::Tau& tau, float default_value) const {
511  get(dnn::n_matched_muons) = n_muons;
512  get(dnn::muon_pt) = best_matched_muon != nullptr ? best_matched_muon->p4().pt() : default_value;
513  get(dnn::muon_dEta) = best_matched_muon != nullptr ? dEta(best_matched_muon->p4(), tau.p4()) : default_value;
514  get(dnn::muon_dPhi) = best_matched_muon != nullptr ? dPhi(best_matched_muon->p4(), tau.p4()) : default_value;
515  get(dnn::muon_n_matches_DT_1) = n_matches.at(MuonSubdetId::DT).at(0);
516  get(dnn::muon_n_matches_DT_2) = n_matches.at(MuonSubdetId::DT).at(1);
517  get(dnn::muon_n_matches_DT_3) = n_matches.at(MuonSubdetId::DT).at(2);
518  get(dnn::muon_n_matches_DT_4) = n_matches.at(MuonSubdetId::DT).at(3);
519  get(dnn::muon_n_matches_CSC_1) = n_matches.at(MuonSubdetId::CSC).at(0);
520  get(dnn::muon_n_matches_CSC_2) = n_matches.at(MuonSubdetId::CSC).at(1);
521  get(dnn::muon_n_matches_CSC_3) = n_matches.at(MuonSubdetId::CSC).at(2);
522  get(dnn::muon_n_matches_CSC_4) = n_matches.at(MuonSubdetId::CSC).at(3);
523  get(dnn::muon_n_hits_DT_2) = n_hits.at(MuonSubdetId::DT).at(1);
524  get(dnn::muon_n_hits_DT_3) = n_hits.at(MuonSubdetId::DT).at(2);
525  get(dnn::muon_n_hits_DT_4) = n_hits.at(MuonSubdetId::DT).at(3);
526  get(dnn::muon_n_hits_CSC_2) = n_hits.at(MuonSubdetId::CSC).at(1);
527  get(dnn::muon_n_hits_CSC_3) = n_hits.at(MuonSubdetId::CSC).at(2);
528  get(dnn::muon_n_hits_CSC_4) = n_hits.at(MuonSubdetId::CSC).at(3);
529  get(dnn::muon_n_hits_RPC_2) = n_hits.at(MuonSubdetId::RPC).at(1);
530  get(dnn::muon_n_hits_RPC_3) = n_hits.at(MuonSubdetId::RPC).at(2);
531  get(dnn::muon_n_hits_RPC_4) = n_hits.at(MuonSubdetId::RPC).at(3);
532  get(dnn::muon_n_stations_with_matches_03) = countMuonStationsWithMatches(0, 3);
533  get(dnn::muon_n_stations_with_hits_23) = countMuonStationsWithHits(2, 3);
534  }
535 
536  private:
537  unsigned countMuonStationsWithMatches(size_t first_station, size_t last_station) const {
538  static const std::map<int, std::vector<bool>> masks = {
539  {MuonSubdetId::DT, {false, false, false, false}},
540  {MuonSubdetId::CSC, {true, false, false, false}},
541  {MuonSubdetId::RPC, {false, false, false, false}},
542  };
543  unsigned cnt = 0;
544  for (unsigned n = first_station; n <= last_station; ++n) {
545  for (const auto& match : n_matches) {
546  if (!masks.at(match.first).at(n) && match.second.at(n) > 0)
547  ++cnt;
548  }
549  }
550  return cnt;
551  }
552 
553  unsigned countMuonStationsWithHits(size_t first_station, size_t last_station) const {
554  static const std::map<int, std::vector<bool>> masks = {
555  {MuonSubdetId::DT, {false, false, false, false}},
556  {MuonSubdetId::CSC, {false, false, false, false}},
557  {MuonSubdetId::RPC, {false, false, false, false}},
558  };
559 
560  unsigned cnt = 0;
561  for (unsigned n = first_station; n <= last_station; ++n) {
562  for (const auto& hit : n_hits) {
563  if (!masks.at(hit.first).at(n) && hit.second.at(n) > 0)
564  ++cnt;
565  }
566  }
567  return cnt;
568  }
569  };
570 
571  struct MuonHitMatchV2 {
572  static constexpr size_t n_muon_stations = 4;
573  static constexpr int first_station_id = 1;
574  static constexpr int last_station_id = first_station_id + n_muon_stations - 1;
575  using CountArray = std::array<unsigned, n_muon_stations>;
576  using CountMap = std::map<int, CountArray>;
577 
578  const std::vector<int>& consideredSubdets() {
579  static const std::vector<int> subdets = {MuonSubdetId::DT, MuonSubdetId::CSC, MuonSubdetId::RPC};
580  return subdets;
581  }
582 
583  const std::string& subdetName(int subdet) {
584  static const std::map<int, std::string> subdet_names = {
585  {MuonSubdetId::DT, "DT"}, {MuonSubdetId::CSC, "CSC"}, {MuonSubdetId::RPC, "RPC"}};
586  if (!subdet_names.count(subdet))
587  throw cms::Exception("MuonHitMatch") << "Subdet name for subdet id " << subdet << " not found.";
588  return subdet_names.at(subdet);
589  }
590 
591  size_t getStationIndex(int station, bool throw_exception) const {
592  if (station < first_station_id || station > last_station_id) {
593  if (throw_exception)
594  throw cms::Exception("MuonHitMatch") << "Station id is out of range";
596  }
597  return static_cast<size_t>(station - 1);
598  }
599 
600  MuonHitMatchV2(const pat::Muon& muon) {
601  for (int subdet : consideredSubdets()) {
602  n_matches[subdet].fill(0);
603  n_hits[subdet].fill(0);
604  }
605 
606  countMatches(muon, n_matches);
607  countHits(muon, n_hits);
608  }
609 
610  void countMatches(const pat::Muon& muon, CountMap& n_matches) {
611  for (const auto& segment : muon.matches()) {
612  if (segment.segmentMatches.empty() && segment.rpcMatches.empty())
613  continue;
614  if (n_matches.count(segment.detector())) {
615  const size_t station_index = getStationIndex(segment.station(), true);
616  ++n_matches.at(segment.detector()).at(station_index);
617  }
618  }
619  }
620 
621  void countHits(const pat::Muon& muon, CountMap& n_hits) {
622  if (muon.outerTrack().isNonnull()) {
623  const auto& hit_pattern = muon.outerTrack()->hitPattern();
624  for (int hit_index = 0; hit_index < hit_pattern.numberOfAllHits(reco::HitPattern::TRACK_HITS); ++hit_index) {
625  auto hit_id = hit_pattern.getHitPattern(reco::HitPattern::TRACK_HITS, hit_index);
626  if (hit_id == 0)
627  break;
628  if (hit_pattern.muonHitFilter(hit_id) && (hit_pattern.getHitType(hit_id) == TrackingRecHit::valid ||
629  hit_pattern.getHitType(hit_id) == TrackingRecHit::bad)) {
630  const size_t station_index = getStationIndex(hit_pattern.getMuonStation(hit_id), false);
631  if (station_index < n_muon_stations) {
632  CountArray* muon_n_hits = nullptr;
633  if (hit_pattern.muonDTHitFilter(hit_id))
634  muon_n_hits = &n_hits.at(MuonSubdetId::DT);
635  else if (hit_pattern.muonCSCHitFilter(hit_id))
636  muon_n_hits = &n_hits.at(MuonSubdetId::CSC);
637  else if (hit_pattern.muonRPCHitFilter(hit_id))
638  muon_n_hits = &n_hits.at(MuonSubdetId::RPC);
639 
640  if (muon_n_hits)
641  ++muon_n_hits->at(station_index);
642  }
643  }
644  }
645  }
646  }
647 
648  unsigned nMatches(int subdet, int station) const {
649  if (!n_matches.count(subdet))
650  throw cms::Exception("MuonHitMatch") << "Subdet " << subdet << " not found.";
651  const size_t station_index = getStationIndex(station, true);
652  return n_matches.at(subdet).at(station_index);
653  }
654 
655  unsigned nHits(int subdet, int station) const {
656  if (!n_hits.count(subdet))
657  throw cms::Exception("MuonHitMatch") << "Subdet " << subdet << " not found.";
658  const size_t station_index = getStationIndex(station, true);
659  return n_hits.at(subdet).at(station_index);
660  }
661 
662  unsigned countMuonStationsWithMatches(int first_station, int last_station) const {
663  static const std::map<int, std::vector<bool>> masks = {
664  {MuonSubdetId::DT, {false, false, false, false}},
665  {MuonSubdetId::CSC, {true, false, false, false}},
666  {MuonSubdetId::RPC, {false, false, false, false}},
667  };
668  const size_t first_station_index = getStationIndex(first_station, true);
669  const size_t last_station_index = getStationIndex(last_station, true);
670  unsigned cnt = 0;
671  for (size_t n = first_station_index; n <= last_station_index; ++n) {
672  for (const auto& match : n_matches) {
673  if (!masks.at(match.first).at(n) && match.second.at(n) > 0)
674  ++cnt;
675  }
676  }
677  return cnt;
678  }
679 
680  unsigned countMuonStationsWithHits(int first_station, int last_station) const {
681  static const std::map<int, std::vector<bool>> masks = {
682  {MuonSubdetId::DT, {false, false, false, false}},
683  {MuonSubdetId::CSC, {false, false, false, false}},
684  {MuonSubdetId::RPC, {false, false, false, false}},
685  };
686 
687  const size_t first_station_index = getStationIndex(first_station, true);
688  const size_t last_station_index = getStationIndex(last_station, true);
689  unsigned cnt = 0;
690  for (size_t n = first_station_index; n <= last_station_index; ++n) {
691  for (const auto& hit : n_hits) {
692  if (!masks.at(hit.first).at(n) && hit.second.at(n) > 0)
693  ++cnt;
694  }
695  }
696  return cnt;
697  }
698 
699  private:
700  CountMap n_matches, n_hits;
701  };
702 
703  enum class CellObjectType {
704  PfCand_electron,
705  PfCand_muon,
706  PfCand_chargedHadron,
707  PfCand_neutralHadron,
708  PfCand_gamma,
709  Electron,
710  Muon,
711  Other
712  };
713 
714  template <typename Object>
715  CellObjectType GetCellObjectType(const Object&);
716  template <>
717  CellObjectType GetCellObjectType(const pat::Electron&) {
719  }
720  template <>
721  CellObjectType GetCellObjectType(const pat::Muon&) {
722  return CellObjectType::Muon;
723  }
724 
725  template <>
726  CellObjectType GetCellObjectType(const pat::PackedCandidate& cand) {
727  static const std::map<int, CellObjectType> obj_types = {{11, CellObjectType::PfCand_electron},
728  {13, CellObjectType::PfCand_muon},
729  {22, CellObjectType::PfCand_gamma},
730  {130, CellObjectType::PfCand_neutralHadron},
731  {211, CellObjectType::PfCand_chargedHadron}};
732 
733  auto iter = obj_types.find(std::abs(cand.pdgId()));
734  if (iter == obj_types.end())
735  return CellObjectType::Other;
736  return iter->second;
737  }
738 
739  using Cell = std::map<CellObjectType, size_t>;
740  struct CellIndex {
741  int eta, phi;
742 
743  bool operator<(const CellIndex& other) const {
744  if (eta != other.eta)
745  return eta < other.eta;
746  return phi < other.phi;
747  }
748  };
749 
750  class CellGrid {
751  public:
752  using Map = std::map<CellIndex, Cell>;
753  using const_iterator = Map::const_iterator;
754 
755  CellGrid(unsigned n_cells_eta, unsigned n_cells_phi, double cell_size_eta, double cell_size_phi)
756  : nCellsEta(n_cells_eta),
757  nCellsPhi(n_cells_phi),
758  nTotal(nCellsEta * nCellsPhi),
759  cellSizeEta(cell_size_eta),
760  cellSizePhi(cell_size_phi) {
761  if (nCellsEta % 2 != 1 || nCellsEta < 1)
762  throw cms::Exception("DeepTauId") << "Invalid number of eta cells.";
763  if (nCellsPhi % 2 != 1 || nCellsPhi < 1)
764  throw cms::Exception("DeepTauId") << "Invalid number of phi cells.";
765  if (cellSizeEta <= 0 || cellSizePhi <= 0)
766  throw cms::Exception("DeepTauId") << "Invalid cell size.";
767  }
768 
769  int maxEtaIndex() const { return static_cast<int>((nCellsEta - 1) / 2); }
770  int maxPhiIndex() const { return static_cast<int>((nCellsPhi - 1) / 2); }
771  double maxDeltaEta() const { return cellSizeEta * (0.5 + maxEtaIndex()); }
772  double maxDeltaPhi() const { return cellSizePhi * (0.5 + maxPhiIndex()); }
773  int getEtaTensorIndex(const CellIndex& cellIndex) const { return cellIndex.eta + maxEtaIndex(); }
774  int getPhiTensorIndex(const CellIndex& cellIndex) const { return cellIndex.phi + maxPhiIndex(); }
775 
776  bool tryGetCellIndex(double deltaEta, double deltaPhi, CellIndex& cellIndex) const {
777  static auto getCellIndex = [](double x, double maxX, double size, int& index) {
778  const double absX = std::abs(x);
779  if (absX > maxX)
780  return false;
781  const double absIndex = std::floor(std::abs(absX / size - 0.5));
782  index = static_cast<int>(std::copysign(absIndex, x));
783  return true;
784  };
785 
786  return getCellIndex(deltaEta, maxDeltaEta(), cellSizeEta, cellIndex.eta) &&
787  getCellIndex(deltaPhi, maxDeltaPhi(), cellSizePhi, cellIndex.phi);
788  }
789 
790  size_t num_valid_cells() const { return cells.size(); }
791  Cell& operator[](const CellIndex& cellIndex) { return cells[cellIndex]; }
792  const Cell& at(const CellIndex& cellIndex) const { return cells.at(cellIndex); }
793  size_t count(const CellIndex& cellIndex) const { return cells.count(cellIndex); }
794  const_iterator find(const CellIndex& cellIndex) const { return cells.find(cellIndex); }
795  const_iterator begin() const { return cells.begin(); }
796  const_iterator end() const { return cells.end(); }
797 
798  public:
799  const unsigned nCellsEta, nCellsPhi, nTotal;
800  const double cellSizeEta, cellSizePhi;
801 
802  private:
803  std::map<CellIndex, Cell> cells;
804  };
805 
806 } // anonymous namespace
807 
809 public:
810  static constexpr float default_value = -999.;
811 
812  static const OutputCollection& GetOutputs() {
813  static constexpr size_t e_index = 0, mu_index = 1, tau_index = 2, jet_index = 3;
814  static const OutputCollection outputs_ = {
815  {"VSe", Output({tau_index}, {e_index, tau_index})},
816  {"VSmu", Output({tau_index}, {mu_index, tau_index})},
817  {"VSjet", Output({tau_index}, {jet_index, tau_index})},
818  };
819  return outputs_;
820  }
821 
824  desc.add<edm::InputTag>("electrons", edm::InputTag("slimmedElectrons"));
825  desc.add<edm::InputTag>("muons", edm::InputTag("slimmedMuons"));
826  desc.add<edm::InputTag>("taus", edm::InputTag("slimmedTaus"));
827  desc.add<edm::InputTag>("pfcands", edm::InputTag("packedPFCandidates"));
828  desc.add<edm::InputTag>("vertices", edm::InputTag("offlineSlimmedPrimaryVertices"));
829  desc.add<edm::InputTag>("rho", edm::InputTag("fixedGridRhoAll"));
830  desc.add<std::vector<std::string>>("graph_file",
831  {"RecoTauTag/TrainingFiles/data/DeepTauId/deepTau_2017v2p6_e6.pb"});
832  desc.add<bool>("mem_mapped", false);
833  desc.add<unsigned>("version", 2);
834  desc.add<int>("debug_level", 0);
835  desc.add<bool>("disable_dxy_pca", false);
836 
837  desc.add<std::vector<std::string>>("VSeWP");
838  desc.add<std::vector<std::string>>("VSmuWP");
839  desc.add<std::vector<std::string>>("VSjetWP");
840  descriptions.add("DeepTau", desc);
841  }
842 
843 public:
846  electrons_token_(consumes<ElectronCollection>(cfg.getParameter<edm::InputTag>("electrons"))),
847  muons_token_(consumes<MuonCollection>(cfg.getParameter<edm::InputTag>("muons"))),
848  rho_token_(consumes<double>(cfg.getParameter<edm::InputTag>("rho"))),
849  version(cfg.getParameter<unsigned>("version")),
850  debug_level(cfg.getParameter<int>("debug_level")),
851  disable_dxy_pca_(cfg.getParameter<bool>("disable_dxy_pca")) {
852  if (version == 1) {
853  input_layer_ = cache_->getGraph().node(0).name();
854  output_layer_ = cache_->getGraph().node(cache_->getGraph().node_size() - 1).name();
855  const auto& shape = cache_->getGraph().node(0).attr().at("shape").shape();
856  if (shape.dim(1).size() != dnn_inputs_2017v1::NumberOfInputs)
857  throw cms::Exception("DeepTauId")
858  << "number of inputs does not match the expected inputs for the given version";
859  } else if (version == 2) {
860  tauBlockTensor_ = std::make_unique<tensorflow::Tensor>(
861  tensorflow::DT_FLOAT, tensorflow::TensorShape{1, dnn_inputs_2017_v2::TauBlockInputs::NumberOfInputs});
862  for (size_t n = 0; n < 2; ++n) {
863  const bool is_inner = n == 0;
864  const auto n_cells =
865  is_inner ? dnn_inputs_2017_v2::number_of_inner_cell : dnn_inputs_2017_v2::number_of_outer_cell;
866  eGammaTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
867  tensorflow::DT_FLOAT,
868  tensorflow::TensorShape{1, 1, 1, dnn_inputs_2017_v2::EgammaBlockInputs::NumberOfInputs});
869  muonTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
870  tensorflow::DT_FLOAT,
871  tensorflow::TensorShape{1, 1, 1, dnn_inputs_2017_v2::MuonBlockInputs::NumberOfInputs});
872  hadronsTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
873  tensorflow::DT_FLOAT,
874  tensorflow::TensorShape{1, 1, 1, dnn_inputs_2017_v2::HadronBlockInputs::NumberOfInputs});
875  convTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
876  tensorflow::DT_FLOAT,
877  tensorflow::TensorShape{1, n_cells, n_cells, dnn_inputs_2017_v2::number_of_conv_features});
878  zeroOutputTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
879  tensorflow::DT_FLOAT, tensorflow::TensorShape{1, 1, 1, dnn_inputs_2017_v2::number_of_conv_features});
880 
881  eGammaTensor_[is_inner]->flat<float>().setZero();
882  muonTensor_[is_inner]->flat<float>().setZero();
883  hadronsTensor_[is_inner]->flat<float>().setZero();
884 
885  setCellConvFeatures(*zeroOutputTensor_[is_inner], getPartialPredictions(is_inner), 0, 0, 0);
886  }
887  } else {
888  throw cms::Exception("DeepTauId") << "version " << version << " is not supported.";
889  }
890  }
891 
892  static std::unique_ptr<deep_tau::DeepTauCache> initializeGlobalCache(const edm::ParameterSet& cfg) {
893  return DeepTauBase::initializeGlobalCache(cfg);
894  }
895 
896  static void globalEndJob(const deep_tau::DeepTauCache* cache_) { return DeepTauBase::globalEndJob(cache_); }
897 
898 private:
899  static constexpr float pi = M_PI;
900 
901  template <typename T>
902  static float getValue(T value) {
903  return std::isnormal(value) ? static_cast<float>(value) : 0.f;
904  }
905 
906  template <typename T>
907  static float getValueLinear(T value, float min_value, float max_value, bool positive) {
908  const float fixed_value = getValue(value);
909  const float clamped_value = std::clamp(fixed_value, min_value, max_value);
910  float transformed_value = (clamped_value - min_value) / (max_value - min_value);
911  if (!positive)
912  transformed_value = transformed_value * 2 - 1;
913  return transformed_value;
914  }
915 
916  template <typename T>
917  static float getValueNorm(T value, float mean, float sigma, float n_sigmas_max = 5) {
918  const float fixed_value = getValue(value);
919  const float norm_value = (fixed_value - mean) / sigma;
920  return std::clamp(norm_value, -n_sigmas_max, n_sigmas_max);
921  }
922 
924  float& cc_ele_energy,
925  float& cc_gamma_energy,
926  int& cc_n_gamma) {
927  cc_ele_energy = cc_gamma_energy = 0;
928  cc_n_gamma = 0;
929  const auto& superCluster = ele.superCluster();
930  if (superCluster.isNonnull() && superCluster.isAvailable() && superCluster->clusters().isNonnull() &&
931  superCluster->clusters().isAvailable()) {
932  for (auto iter = superCluster->clustersBegin(); iter != superCluster->clustersEnd(); ++iter) {
933  const float energy = static_cast<float>((*iter)->energy());
934  if (iter == superCluster->clustersBegin())
935  cc_ele_energy += energy;
936  else {
937  cc_gamma_energy += energy;
938  ++cc_n_gamma;
939  }
940  }
941  return true;
942  } else
943  return false;
944  }
945 
946  inline void checkInputs(
947  const tensorflow::Tensor& inputs, const char* block_name, int n_inputs, int n_eta = 1, int n_phi = 1) const {
948  if (debug_level >= 1) {
949  for (int eta = 0; eta < n_eta; ++eta) {
950  for (int phi = 0; phi < n_phi; phi++) {
951  for (int k = 0; k < n_inputs; ++k) {
952  const float input =
953  n_eta == 1 && n_phi == 1 ? inputs.matrix<float>()(0, k) : inputs.tensor<float, 4>()(0, eta, phi, k);
954  if (edm::isNotFinite(input))
955  throw cms::Exception("DeepTauId")
956  << "in the " << block_name << ", input is not finite, i.e. infinite or NaN, for eta_index = " << n_eta
957  << ", phi_index = " << n_phi << ", input_index = " << k;
958  if (debug_level >= 2)
959  std::cout << block_name << "," << eta << "," << phi << "," << k << "," << std::setprecision(5)
960  << std::fixed << input << '\n';
961  }
962  }
963  }
964  }
965  }
966 
967 private:
968  tensorflow::Tensor getPredictions(edm::Event& event,
969  const edm::EventSetup& es,
970  edm::Handle<TauCollection> taus) override {
972  event.getByToken(electrons_token_, electrons);
973 
975  event.getByToken(muons_token_, muons);
976 
978  event.getByToken(pfcandToken_, pfCands);
979 
981  event.getByToken(vtxToken_, vertices);
982 
984  event.getByToken(rho_token_, rho);
985 
986  tensorflow::Tensor predictions(tensorflow::DT_FLOAT, {static_cast<int>(taus->size()), deep_tau::NumberOfOutputs});
987  for (size_t tau_index = 0; tau_index < taus->size(); ++tau_index) {
988  std::vector<tensorflow::Tensor> pred_vector;
989  if (version == 1)
990  getPredictionsV1(taus->at(tau_index), *electrons, *muons, pred_vector);
991  else if (version == 2)
992  getPredictionsV2(taus->at(tau_index), *electrons, *muons, *pfCands, vertices->at(0), *rho, pred_vector);
993  else
994  throw cms::Exception("DeepTauId") << "version " << version << " is not supported.";
995  for (int k = 0; k < deep_tau::NumberOfOutputs; ++k) {
996  const float pred = pred_vector[0].flat<float>()(k);
997  if (!(pred >= 0 && pred <= 1))
998  throw cms::Exception("DeepTauId")
999  << "invalid prediction = " << pred << " for tau_index = " << tau_index << ", pred_index = " << k;
1000  predictions.matrix<float>()(tau_index, k) = pred;
1001  }
1002  }
1003  return predictions;
1004  }
1005 
1008  const pat::MuonCollection& muons,
1009  std::vector<tensorflow::Tensor>& pred_vector) {
1010  const tensorflow::Tensor& inputs = createInputsV1<dnn_inputs_2017v1>(tau, electrons, muons);
1011  tensorflow::run(&(cache_->getSession()), {{input_layer_, inputs}}, {output_layer_}, &pred_vector);
1012  }
1013 
1016  const pat::MuonCollection& muons,
1017  const pat::PackedCandidateCollection& pfCands,
1018  const reco::Vertex& pv,
1019  double rho,
1020  std::vector<tensorflow::Tensor>& pred_vector) {
1021  CellGrid inner_grid(dnn_inputs_2017_v2::number_of_inner_cell, dnn_inputs_2017_v2::number_of_inner_cell, 0.02, 0.02);
1022  CellGrid outer_grid(dnn_inputs_2017_v2::number_of_outer_cell, dnn_inputs_2017_v2::number_of_outer_cell, 0.05, 0.05);
1023  fillGrids(tau, electrons, inner_grid, outer_grid);
1024  fillGrids(tau, muons, inner_grid, outer_grid);
1025  fillGrids(tau, pfCands, inner_grid, outer_grid);
1026 
1028  createConvFeatures(tau, pv, rho, electrons, muons, pfCands, inner_grid, true);
1029  createConvFeatures(tau, pv, rho, electrons, muons, pfCands, outer_grid, false);
1030 
1031  tensorflow::run(&(cache_->getSession("core")),
1032  {{"input_tau", *tauBlockTensor_},
1033  {"input_inner", *convTensor_.at(true)},
1034  {"input_outer", *convTensor_.at(false)}},
1035  {"main_output/Softmax"},
1036  &pred_vector);
1037  }
1038 
1039  template <typename Collection>
1040  void fillGrids(const TauType& tau, const Collection& objects, CellGrid& inner_grid, CellGrid& outer_grid) {
1041  static constexpr double outer_dR2 = 0.25; //0.5^2
1042  const double inner_radius = getInnerSignalConeRadius(tau.polarP4().pt());
1043  const double inner_dR2 = std::pow(inner_radius, 2);
1044 
1045  const auto addObject = [&](size_t n, double deta, double dphi, CellGrid& grid) {
1046  const auto& obj = objects.at(n);
1047  const CellObjectType obj_type = GetCellObjectType(obj);
1048  if (obj_type == CellObjectType::Other)
1049  return;
1050  CellIndex cell_index;
1051  if (grid.tryGetCellIndex(deta, dphi, cell_index)) {
1052  Cell& cell = grid[cell_index];
1053  auto iter = cell.find(obj_type);
1054  if (iter != cell.end()) {
1055  const auto& prev_obj = objects.at(iter->second);
1056  if (obj.polarP4().pt() > prev_obj.polarP4().pt())
1057  iter->second = n;
1058  } else {
1059  cell[obj_type] = n;
1060  }
1061  }
1062  };
1063 
1064  for (size_t n = 0; n < objects.size(); ++n) {
1065  const auto& obj = objects.at(n);
1066  const double deta = obj.polarP4().eta() - tau.polarP4().eta();
1067  const double dphi = reco::deltaPhi(obj.polarP4().phi(), tau.polarP4().phi());
1068  const double dR2 = std::pow(deta, 2) + std::pow(dphi, 2);
1069  if (dR2 < inner_dR2)
1070  addObject(n, deta, dphi, inner_grid);
1071  if (dR2 < outer_dR2)
1072  addObject(n, deta, dphi, outer_grid);
1073  }
1074  }
1075 
1076  tensorflow::Tensor getPartialPredictions(bool is_inner) {
1077  std::vector<tensorflow::Tensor> pred_vector;
1078  if (is_inner) {
1079  tensorflow::run(&(cache_->getSession("inner")),
1080  {
1081  {"input_inner_egamma", *eGammaTensor_.at(is_inner)},
1082  {"input_inner_muon", *muonTensor_.at(is_inner)},
1083  {"input_inner_hadrons", *hadronsTensor_.at(is_inner)},
1084  },
1085  {"inner_all_dropout_4/Identity"},
1086  &pred_vector);
1087  } else {
1088  tensorflow::run(&(cache_->getSession("outer")),
1089  {
1090  {"input_outer_egamma", *eGammaTensor_.at(is_inner)},
1091  {"input_outer_muon", *muonTensor_.at(is_inner)},
1092  {"input_outer_hadrons", *hadronsTensor_.at(is_inner)},
1093  },
1094  {"outer_all_dropout_4/Identity"},
1095  &pred_vector);
1096  }
1097  return pred_vector.at(0);
1098  }
1099 
1101  const reco::Vertex& pv,
1102  double rho,
1104  const pat::MuonCollection& muons,
1105  const pat::PackedCandidateCollection& pfCands,
1106  const CellGrid& grid,
1107  bool is_inner) {
1108  tensorflow::Tensor& convTensor = *convTensor_.at(is_inner);
1109 
1110  eGammaTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1111  tensorflow::DT_FLOAT,
1112  tensorflow::TensorShape{
1113  (long long int)grid.num_valid_cells(), 1, 1, dnn_inputs_2017_v2::EgammaBlockInputs::NumberOfInputs});
1114  muonTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1115  tensorflow::DT_FLOAT,
1116  tensorflow::TensorShape{
1117  (long long int)grid.num_valid_cells(), 1, 1, dnn_inputs_2017_v2::MuonBlockInputs::NumberOfInputs});
1118  hadronsTensor_[is_inner] = std::make_unique<tensorflow::Tensor>(
1119  tensorflow::DT_FLOAT,
1120  tensorflow::TensorShape{
1121  (long long int)grid.num_valid_cells(), 1, 1, dnn_inputs_2017_v2::HadronBlockInputs::NumberOfInputs});
1122 
1123  eGammaTensor_[is_inner]->flat<float>().setZero();
1124  muonTensor_[is_inner]->flat<float>().setZero();
1125  hadronsTensor_[is_inner]->flat<float>().setZero();
1126 
1127  unsigned idx = 0;
1128  for (int eta = -grid.maxEtaIndex(); eta <= grid.maxEtaIndex(); ++eta) {
1129  for (int phi = -grid.maxPhiIndex(); phi <= grid.maxPhiIndex(); ++phi) {
1130  const CellIndex cell_index{eta, phi};
1131  const auto cell_iter = grid.find(cell_index);
1132  if (cell_iter != grid.end()) {
1133  const Cell& cell = cell_iter->second;
1134  createEgammaBlockInputs(idx, tau, pv, rho, electrons, pfCands, cell, is_inner);
1135  createMuonBlockInputs(idx, tau, pv, rho, muons, pfCands, cell, is_inner);
1136  createHadronsBlockInputs(idx, tau, pv, rho, pfCands, cell, is_inner);
1137  idx += 1;
1138  }
1139  }
1140  }
1141 
1142  const auto predTensor = getPartialPredictions(is_inner);
1143 
1144  idx = 0;
1145  for (int eta = -grid.maxEtaIndex(); eta <= grid.maxEtaIndex(); ++eta) {
1146  for (int phi = -grid.maxPhiIndex(); phi <= grid.maxPhiIndex(); ++phi) {
1147  const CellIndex cell_index{eta, phi};
1148  const int eta_index = grid.getEtaTensorIndex(cell_index);
1149  const int phi_index = grid.getPhiTensorIndex(cell_index);
1150 
1151  const auto cell_iter = grid.find(cell_index);
1152  if (cell_iter != grid.end()) {
1153  setCellConvFeatures(convTensor, predTensor, idx, eta_index, phi_index);
1154  idx += 1;
1155  } else {
1156  setCellConvFeatures(convTensor, *zeroOutputTensor_[is_inner], 0, eta_index, phi_index);
1157  }
1158  }
1159  }
1160  }
1161 
1162  void setCellConvFeatures(tensorflow::Tensor& convTensor,
1163  const tensorflow::Tensor& features,
1164  unsigned batch_idx,
1165  int eta_index,
1166  int phi_index) {
1167  for (int n = 0; n < dnn_inputs_2017_v2::number_of_conv_features; ++n)
1168  convTensor.tensor<float, 4>()(0, eta_index, phi_index, n) = features.tensor<float, 4>()(batch_idx, 0, 0, n);
1169  }
1170 
1171  void createTauBlockInputs(const TauType& tau, const reco::Vertex& pv, double rho) {
1172  namespace dnn = dnn_inputs_2017_v2::TauBlockInputs;
1173 
1174  tensorflow::Tensor& inputs = *tauBlockTensor_;
1175  inputs.flat<float>().setZero();
1176 
1177  const auto& get = [&](int var_index) -> float& { return inputs.matrix<float>()(0, var_index); };
1178 
1179  auto leadChargedHadrCand = dynamic_cast<const pat::PackedCandidate*>(tau.leadChargedHadrCand().get());
1180 
1181  get(dnn::rho) = getValueNorm(rho, 21.49f, 9.713f);
1182  get(dnn::tau_pt) = getValueLinear(tau.polarP4().pt(), 20.f, 1000.f, true);
1183  get(dnn::tau_eta) = getValueLinear(tau.polarP4().eta(), -2.3f, 2.3f, false);
1184  get(dnn::tau_phi) = getValueLinear(tau.polarP4().phi(), -pi, pi, false);
1185  get(dnn::tau_mass) = getValueNorm(tau.polarP4().mass(), 0.6669f, 0.6553f);
1186  get(dnn::tau_E_over_pt) = getValueLinear(tau.p4().energy() / tau.p4().pt(), 1.f, 5.2f, true);
1187  get(dnn::tau_charge) = getValue(tau.charge());
1188  get(dnn::tau_n_charged_prongs) = getValueLinear(tau.decayMode() / 5 + 1, 1, 3, true);
1189  get(dnn::tau_n_neutral_prongs) = getValueLinear(tau.decayMode() % 5, 0, 2, true);
1190  get(dnn::chargedIsoPtSum) = getValueNorm(tau.tauID("chargedIsoPtSum"), 47.78f, 123.5f);
1191  get(dnn::chargedIsoPtSumdR03_over_dR05) = getValue(tau.tauID("chargedIsoPtSumdR03") / tau.tauID("chargedIsoPtSum"));
1192  get(dnn::footprintCorrection) = getValueNorm(tau.tauID("footprintCorrectiondR03"), 9.029f, 26.42f);
1193  get(dnn::neutralIsoPtSum) = getValueNorm(tau.tauID("neutralIsoPtSum"), 57.59f, 155.3f);
1194  get(dnn::neutralIsoPtSumWeight_over_neutralIsoPtSum) =
1195  getValue(tau.tauID("neutralIsoPtSumWeight") / tau.tauID("neutralIsoPtSum"));
1196  get(dnn::neutralIsoPtSumWeightdR03_over_neutralIsoPtSum) =
1197  getValue(tau.tauID("neutralIsoPtSumWeightdR03") / tau.tauID("neutralIsoPtSum"));
1198  get(dnn::neutralIsoPtSumdR03_over_dR05) = getValue(tau.tauID("neutralIsoPtSumdR03") / tau.tauID("neutralIsoPtSum"));
1199  get(dnn::photonPtSumOutsideSignalCone) =
1200  getValueNorm(tau.tauID("photonPtSumOutsideSignalConedR03"), 1.731f, 6.846f);
1201  get(dnn::puCorrPtSum) = getValueNorm(tau.tauID("puCorrPtSum"), 22.38f, 16.34f);
1202  // The global PCA coordinates were used as inputs during the NN training, but it was decided to disable
1203  // them for the inference, because modeling of dxy_PCA in MC poorly describes the data, and x and y coordinates
1204  // in data results outside of the expected 5 std. dev. input validity range. On the other hand,
1205  // these coordinates are strongly era-dependent. Kept as comment to document what NN expects.
1206  if (!disable_dxy_pca_) {
1207  get(dnn::tau_dxy_pca_x) = getValueNorm(tau.dxy_PCA().x(), -0.0241f, 0.0074f);
1208  get(dnn::tau_dxy_pca_y) = getValueNorm(tau.dxy_PCA().y(), 0.0675f, 0.0128f);
1209  get(dnn::tau_dxy_pca_z) = getValueNorm(tau.dxy_PCA().z(), 0.7973f, 3.456f);
1210  } else {
1211  get(dnn::tau_dxy_pca_x) = 0;
1212  get(dnn::tau_dxy_pca_y) = 0;
1213  get(dnn::tau_dxy_pca_z) = 0;
1214  }
1215 
1216  const bool tau_dxy_valid =
1217  std::isnormal(tau.dxy()) && tau.dxy() > -10 && std::isnormal(tau.dxy_error()) && tau.dxy_error() > 0;
1218  if (tau_dxy_valid) {
1219  get(dnn::tau_dxy_valid) = tau_dxy_valid;
1220  get(dnn::tau_dxy) = getValueNorm(tau.dxy(), 0.0018f, 0.0085f);
1221  get(dnn::tau_dxy_sig) = getValueNorm(std::abs(tau.dxy()) / tau.dxy_error(), 2.26f, 4.191f);
1222  }
1223  const bool tau_ip3d_valid =
1224  std::isnormal(tau.ip3d()) && tau.ip3d() > -10 && std::isnormal(tau.ip3d_error()) && tau.ip3d_error() > 0;
1225  if (tau_ip3d_valid) {
1226  get(dnn::tau_ip3d_valid) = tau_ip3d_valid;
1227  get(dnn::tau_ip3d) = getValueNorm(tau.ip3d(), 0.0026f, 0.0114f);
1228  get(dnn::tau_ip3d_sig) = getValueNorm(std::abs(tau.ip3d()) / tau.ip3d_error(), 2.928f, 4.466f);
1229  }
1230  if (leadChargedHadrCand) {
1231  get(dnn::tau_dz) = getValueNorm(leadChargedHadrCand->dz(), 0.f, 0.0190f);
1232  const bool tau_dz_sig_valid = leadChargedHadrCand->hasTrackDetails() &&
1233  std::isnormal(leadChargedHadrCand->dz()) &&
1234  std::isnormal(leadChargedHadrCand->dzError()) && leadChargedHadrCand->dzError() > 0;
1235  get(dnn::tau_dz_sig_valid) = tau_dz_sig_valid;
1236  const double dzError = leadChargedHadrCand->hasTrackDetails() ? leadChargedHadrCand->dzError() : default_value;
1237  get(dnn::tau_dz_sig) = getValueNorm(std::abs(leadChargedHadrCand->dz()) / dzError, 4.717f, 11.78f);
1238  }
1239  get(dnn::tau_flightLength_x) = getValueNorm(tau.flightLength().x(), -0.0003f, 0.7362f);
1240  get(dnn::tau_flightLength_y) = getValueNorm(tau.flightLength().y(), -0.0009f, 0.7354f);
1241  get(dnn::tau_flightLength_z) = getValueNorm(tau.flightLength().z(), -0.0022f, 1.993f);
1242  // get(dnn::tau_flightLength_sig) = getValueNorm(tau.flightLengthSig(), -4.78f, 9.573f);
1243  get(dnn::tau_flightLength_sig) = 0.55756444; //This value is set due to a bug in the training
1244  get(dnn::tau_pt_weighted_deta_strip) =
1245  getValueLinear(reco::tau::pt_weighted_deta_strip(tau, tau.decayMode()), 0, 1, true);
1246 
1247  get(dnn::tau_pt_weighted_dphi_strip) =
1248  getValueLinear(reco::tau::pt_weighted_dphi_strip(tau, tau.decayMode()), 0, 1, true);
1249  get(dnn::tau_pt_weighted_dr_signal) =
1250  getValueNorm(reco::tau::pt_weighted_dr_signal(tau, tau.decayMode()), 0.0052f, 0.01433f);
1251  get(dnn::tau_pt_weighted_dr_iso) = getValueLinear(reco::tau::pt_weighted_dr_iso(tau, tau.decayMode()), 0, 1, true);
1252  get(dnn::tau_leadingTrackNormChi2) = getValueNorm(tau.leadingTrackNormChi2(), 1.538f, 4.401f);
1253  const auto eratio = reco::tau::eratio(tau);
1254  const bool tau_e_ratio_valid = std::isnormal(eratio) && eratio > 0.f;
1255  get(dnn::tau_e_ratio_valid) = tau_e_ratio_valid;
1256  get(dnn::tau_e_ratio) = tau_e_ratio_valid ? getValueLinear(eratio, 0, 1, true) : 0.f;
1257  const double gj_angle_diff = calculateGottfriedJacksonAngleDifference(tau);
1258  const bool tau_gj_angle_diff_valid = (std::isnormal(gj_angle_diff) || gj_angle_diff == 0) && gj_angle_diff >= 0;
1259  get(dnn::tau_gj_angle_diff_valid) = tau_gj_angle_diff_valid;
1260  get(dnn::tau_gj_angle_diff) = tau_gj_angle_diff_valid ? getValueLinear(gj_angle_diff, 0, pi, true) : 0;
1261  get(dnn::tau_n_photons) = getValueNorm(reco::tau::n_photons_total(tau), 2.95f, 3.927f);
1262  get(dnn::tau_emFraction) = getValueLinear(tau.emFraction_MVA(), -1, 1, false);
1263  get(dnn::tau_inside_ecal_crack) = getValue(isInEcalCrack(tau.p4().eta()));
1264  get(dnn::leadChargedCand_etaAtEcalEntrance_minus_tau_eta) =
1265  getValueNorm(tau.etaAtEcalEntranceLeadChargedCand() - tau.p4().eta(), 0.0042f, 0.0323f);
1266 
1267  checkInputs(inputs, "tau_block", dnn::NumberOfInputs);
1268  }
1269 
1271  const TauType& tau,
1272  const reco::Vertex& pv,
1273  double rho,
1275  const pat::PackedCandidateCollection& pfCands,
1276  const Cell& cell_map,
1277  bool is_inner) {
1278  namespace dnn = dnn_inputs_2017_v2::EgammaBlockInputs;
1279 
1280  tensorflow::Tensor& inputs = *eGammaTensor_.at(is_inner);
1281 
1282  const auto& get = [&](int var_index) -> float& { return inputs.tensor<float, 4>()(idx, 0, 0, var_index); };
1283 
1284  const bool valid_index_pf_ele = cell_map.count(CellObjectType::PfCand_electron);
1285  const bool valid_index_pf_gamma = cell_map.count(CellObjectType::PfCand_gamma);
1286  const bool valid_index_ele = cell_map.count(CellObjectType::Electron);
1287 
1288  if (!cell_map.empty()) {
1289  get(dnn::rho) = getValueNorm(rho, 21.49f, 9.713f);
1290  get(dnn::tau_pt) = getValueLinear(tau.polarP4().pt(), 20.f, 1000.f, true);
1291  get(dnn::tau_eta) = getValueLinear(tau.polarP4().eta(), -2.3f, 2.3f, false);
1292  get(dnn::tau_inside_ecal_crack) = getValue(isInEcalCrack(tau.polarP4().eta()));
1293  }
1294  if (valid_index_pf_ele) {
1295  size_t index_pf_ele = cell_map.at(CellObjectType::PfCand_electron);
1296 
1297  get(dnn::pfCand_ele_valid) = valid_index_pf_ele;
1298  get(dnn::pfCand_ele_rel_pt) = getValueNorm(pfCands.at(index_pf_ele).polarP4().pt() / tau.polarP4().pt(),
1299  is_inner ? 0.9792f : 0.304f,
1300  is_inner ? 0.5383f : 1.845f);
1301  get(dnn::pfCand_ele_deta) = getValueLinear(pfCands.at(index_pf_ele).polarP4().eta() - tau.polarP4().eta(),
1302  is_inner ? -0.1f : -0.5f,
1303  is_inner ? 0.1f : 0.5f,
1304  false);
1305  get(dnn::pfCand_ele_dphi) = getValueLinear(dPhi(tau.polarP4(), pfCands.at(index_pf_ele).polarP4()),
1306  is_inner ? -0.1f : -0.5f,
1307  is_inner ? 0.1f : 0.5f,
1308  false);
1309  get(dnn::pfCand_ele_pvAssociationQuality) =
1310  getValueLinear<int>(pfCands.at(index_pf_ele).pvAssociationQuality(), 0, 7, true);
1311  get(dnn::pfCand_ele_puppiWeight) = getValue(pfCands.at(index_pf_ele).puppiWeight());
1312  get(dnn::pfCand_ele_charge) = getValue(pfCands.at(index_pf_ele).charge());
1313  get(dnn::pfCand_ele_lostInnerHits) = getValue<int>(pfCands.at(index_pf_ele).lostInnerHits());
1314  get(dnn::pfCand_ele_numberOfPixelHits) =
1315  getValueLinear(pfCands.at(index_pf_ele).numberOfPixelHits(), 0, 10, true);
1316  get(dnn::pfCand_ele_vertex_dx) =
1317  getValueNorm(pfCands.at(index_pf_ele).vertex().x() - pv.position().x(), 0.f, 0.1221f);
1318  get(dnn::pfCand_ele_vertex_dy) =
1319  getValueNorm(pfCands.at(index_pf_ele).vertex().y() - pv.position().y(), 0.f, 0.1226f);
1320  get(dnn::pfCand_ele_vertex_dz) =
1321  getValueNorm(pfCands.at(index_pf_ele).vertex().z() - pv.position().z(), 0.001f, 1.024f);
1322  get(dnn::pfCand_ele_vertex_dx_tauFL) = getValueNorm(
1323  pfCands.at(index_pf_ele).vertex().x() - pv.position().x() - tau.flightLength().x(), 0.f, 0.3411f);
1324  get(dnn::pfCand_ele_vertex_dy_tauFL) = getValueNorm(
1325  pfCands.at(index_pf_ele).vertex().y() - pv.position().y() - tau.flightLength().y(), 0.0003f, 0.3385f);
1326  get(dnn::pfCand_ele_vertex_dz_tauFL) =
1327  getValueNorm(pfCands.at(index_pf_ele).vertex().z() - pv.position().z() - tau.flightLength().z(), 0.f, 1.307f);
1328 
1329  const bool hasTrackDetails = pfCands.at(index_pf_ele).hasTrackDetails();
1330  if (hasTrackDetails) {
1331  get(dnn::pfCand_ele_hasTrackDetails) = hasTrackDetails;
1332  get(dnn::pfCand_ele_dxy) = getValueNorm(pfCands.at(index_pf_ele).dxy(), 0.f, 0.171f);
1333  get(dnn::pfCand_ele_dxy_sig) =
1334  getValueNorm(std::abs(pfCands.at(index_pf_ele).dxy()) / pfCands.at(index_pf_ele).dxyError(), 1.634f, 6.45f);
1335  get(dnn::pfCand_ele_dz) = getValueNorm(pfCands.at(index_pf_ele).dz(), 0.001f, 1.02f);
1336  get(dnn::pfCand_ele_dz_sig) =
1337  getValueNorm(std::abs(pfCands.at(index_pf_ele).dz()) / pfCands.at(index_pf_ele).dzError(), 24.56f, 210.4f);
1338  get(dnn::pfCand_ele_track_chi2_ndof) =
1339  getValueNorm(pfCands.at(index_pf_ele).pseudoTrack().chi2() / pfCands.at(index_pf_ele).pseudoTrack().ndof(),
1340  2.272f,
1341  8.439f);
1342  get(dnn::pfCand_ele_track_ndof) = getValueNorm(pfCands.at(index_pf_ele).pseudoTrack().ndof(), 15.18f, 3.203f);
1343  }
1344  }
1345  if (valid_index_pf_gamma) {
1346  size_t index_pf_gamma = cell_map.at(CellObjectType::PfCand_gamma);
1347  get(dnn::pfCand_gamma_valid) = valid_index_pf_gamma;
1348  get(dnn::pfCand_gamma_rel_pt) = getValueNorm(pfCands.at(index_pf_gamma).polarP4().pt() / tau.polarP4().pt(),
1349  is_inner ? 0.6048f : 0.02576f,
1350  is_inner ? 1.669f : 0.3833f);
1351  get(dnn::pfCand_gamma_deta) = getValueLinear(pfCands.at(index_pf_gamma).polarP4().eta() - tau.polarP4().eta(),
1352  is_inner ? -0.1f : -0.5f,
1353  is_inner ? 0.1f : 0.5f,
1354  false);
1355  get(dnn::pfCand_gamma_dphi) = getValueLinear(dPhi(tau.polarP4(), pfCands.at(index_pf_gamma).polarP4()),
1356  is_inner ? -0.1f : -0.5f,
1357  is_inner ? 0.1f : 0.5f,
1358  false);
1359  get(dnn::pfCand_gamma_pvAssociationQuality) =
1360  getValueLinear<int>(pfCands.at(index_pf_gamma).pvAssociationQuality(), 0, 7, true);
1361  get(dnn::pfCand_gamma_fromPV) = getValueLinear<int>(pfCands.at(index_pf_gamma).fromPV(), 0, 3, true);
1362  get(dnn::pfCand_gamma_puppiWeight) = getValue(pfCands.at(index_pf_gamma).puppiWeight());
1363  get(dnn::pfCand_gamma_puppiWeightNoLep) = getValue(pfCands.at(index_pf_gamma).puppiWeightNoLep());
1364  get(dnn::pfCand_gamma_lostInnerHits) = getValue<int>(pfCands.at(index_pf_gamma).lostInnerHits());
1365  get(dnn::pfCand_gamma_numberOfPixelHits) =
1366  getValueLinear(pfCands.at(index_pf_gamma).numberOfPixelHits(), 0, 7, true);
1367  get(dnn::pfCand_gamma_vertex_dx) =
1368  getValueNorm(pfCands.at(index_pf_gamma).vertex().x() - pv.position().x(), 0.f, 0.0067f);
1369  get(dnn::pfCand_gamma_vertex_dy) =
1370  getValueNorm(pfCands.at(index_pf_gamma).vertex().y() - pv.position().y(), 0.f, 0.0069f);
1371  get(dnn::pfCand_gamma_vertex_dz) =
1372  getValueNorm(pfCands.at(index_pf_gamma).vertex().z() - pv.position().z(), 0.f, 0.0578f);
1373  get(dnn::pfCand_gamma_vertex_dx_tauFL) = getValueNorm(
1374  pfCands.at(index_pf_gamma).vertex().x() - pv.position().x() - tau.flightLength().x(), 0.001f, 0.9565f);
1375  get(dnn::pfCand_gamma_vertex_dy_tauFL) = getValueNorm(
1376  pfCands.at(index_pf_gamma).vertex().y() - pv.position().y() - tau.flightLength().y(), 0.0008f, 0.9592f);
1377  get(dnn::pfCand_gamma_vertex_dz_tauFL) = getValueNorm(
1378  pfCands.at(index_pf_gamma).vertex().z() - pv.position().z() - tau.flightLength().z(), 0.0038f, 2.154f);
1379 
1380  const bool hasTrackDetails = pfCands.at(index_pf_gamma).hasTrackDetails();
1381  if (hasTrackDetails) {
1382  get(dnn::pfCand_gamma_hasTrackDetails) = hasTrackDetails;
1383  get(dnn::pfCand_gamma_dxy) = getValueNorm(pfCands.at(index_pf_gamma).dxy(), 0.0004f, 0.882f);
1384  get(dnn::pfCand_gamma_dxy_sig) = getValueNorm(
1385  std::abs(pfCands.at(index_pf_gamma).dxy()) / pfCands.at(index_pf_gamma).dxyError(), 4.271f, 63.78f);
1386  get(dnn::pfCand_gamma_dz) = getValueNorm(pfCands.at(index_pf_gamma).dz(), 0.0071f, 5.285f);
1387  get(dnn::pfCand_gamma_dz_sig) = getValueNorm(
1388  std::abs(pfCands.at(index_pf_gamma).dz()) / pfCands.at(index_pf_gamma).dzError(), 162.1f, 622.4f);
1389  get(dnn::pfCand_gamma_track_chi2_ndof) = pfCands.at(index_pf_gamma).pseudoTrack().ndof() > 0
1390  ? getValueNorm(pfCands.at(index_pf_gamma).pseudoTrack().chi2() /
1391  pfCands.at(index_pf_gamma).pseudoTrack().ndof(),
1392  4.268f,
1393  15.47f)
1394  : 0;
1395  get(dnn::pfCand_gamma_track_ndof) =
1396  pfCands.at(index_pf_gamma).pseudoTrack().ndof() > 0
1397  ? getValueNorm(pfCands.at(index_pf_gamma).pseudoTrack().ndof(), 12.25f, 4.774f)
1398  : 0;
1399  }
1400  }
1401  if (valid_index_ele) {
1402  size_t index_ele = cell_map.at(CellObjectType::Electron);
1403 
1404  get(dnn::ele_valid) = valid_index_ele;
1405  get(dnn::ele_rel_pt) = getValueNorm(electrons.at(index_ele).polarP4().pt() / tau.polarP4().pt(),
1406  is_inner ? 1.067f : 0.5111f,
1407  is_inner ? 1.521f : 2.765f);
1408  get(dnn::ele_deta) = getValueLinear(electrons.at(index_ele).polarP4().eta() - tau.polarP4().eta(),
1409  is_inner ? -0.1f : -0.5f,
1410  is_inner ? 0.1f : 0.5f,
1411  false);
1412  get(dnn::ele_dphi) = getValueLinear(dPhi(tau.polarP4(), electrons.at(index_ele).polarP4()),
1413  is_inner ? -0.1f : -0.5f,
1414  is_inner ? 0.1f : 0.5f,
1415  false);
1416 
1417  float cc_ele_energy, cc_gamma_energy;
1418  int cc_n_gamma;
1419  const bool cc_valid =
1420  calculateElectronClusterVarsV2(electrons.at(index_ele), cc_ele_energy, cc_gamma_energy, cc_n_gamma);
1421  if (cc_valid) {
1422  get(dnn::ele_cc_valid) = cc_valid;
1423  get(dnn::ele_cc_ele_rel_energy) =
1424  getValueNorm(cc_ele_energy / electrons.at(index_ele).polarP4().pt(), 1.729f, 1.644f);
1425  get(dnn::ele_cc_gamma_rel_energy) = getValueNorm(cc_gamma_energy / cc_ele_energy, 0.1439f, 0.3284f);
1426  get(dnn::ele_cc_n_gamma) = getValueNorm(cc_n_gamma, 1.794f, 2.079f);
1427  }
1428  get(dnn::ele_rel_trackMomentumAtVtx) = getValueNorm(
1429  electrons.at(index_ele).trackMomentumAtVtx().R() / electrons.at(index_ele).polarP4().pt(), 1.531f, 1.424f);
1430  get(dnn::ele_rel_trackMomentumAtCalo) = getValueNorm(
1431  electrons.at(index_ele).trackMomentumAtCalo().R() / electrons.at(index_ele).polarP4().pt(), 1.531f, 1.424f);
1432  get(dnn::ele_rel_trackMomentumOut) = getValueNorm(
1433  electrons.at(index_ele).trackMomentumOut().R() / electrons.at(index_ele).polarP4().pt(), 0.7735f, 0.935f);
1434  get(dnn::ele_rel_trackMomentumAtEleClus) =
1435  getValueNorm(electrons.at(index_ele).trackMomentumAtEleClus().R() / electrons.at(index_ele).polarP4().pt(),
1436  0.7735f,
1437  0.935f);
1438  get(dnn::ele_rel_trackMomentumAtVtxWithConstraint) = getValueNorm(
1439  electrons.at(index_ele).trackMomentumAtVtxWithConstraint().R() / electrons.at(index_ele).polarP4().pt(),
1440  1.625f,
1441  1.581f);
1442  get(dnn::ele_rel_ecalEnergy) =
1443  getValueNorm(electrons.at(index_ele).ecalEnergy() / electrons.at(index_ele).polarP4().pt(), 1.993f, 1.308f);
1444  get(dnn::ele_ecalEnergy_sig) = getValueNorm(
1445  electrons.at(index_ele).ecalEnergy() / electrons.at(index_ele).ecalEnergyError(), 70.25f, 58.16f);
1446  get(dnn::ele_eSuperClusterOverP) = getValueNorm(electrons.at(index_ele).eSuperClusterOverP(), 2.432f, 15.13f);
1447  get(dnn::ele_eSeedClusterOverP) = getValueNorm(electrons.at(index_ele).eSeedClusterOverP(), 2.034f, 13.96f);
1448  get(dnn::ele_eSeedClusterOverPout) = getValueNorm(electrons.at(index_ele).eSeedClusterOverPout(), 6.64f, 36.8f);
1449  get(dnn::ele_eEleClusterOverPout) = getValueNorm(electrons.at(index_ele).eEleClusterOverPout(), 4.183f, 20.63f);
1450  get(dnn::ele_deltaEtaSuperClusterTrackAtVtx) =
1451  getValueNorm(electrons.at(index_ele).deltaEtaSuperClusterTrackAtVtx(), 0.f, 0.0363f);
1452  get(dnn::ele_deltaEtaSeedClusterTrackAtCalo) =
1453  getValueNorm(electrons.at(index_ele).deltaEtaSeedClusterTrackAtCalo(), -0.0001f, 0.0512f);
1454  get(dnn::ele_deltaEtaEleClusterTrackAtCalo) =
1455  getValueNorm(electrons.at(index_ele).deltaEtaEleClusterTrackAtCalo(), -0.0001f, 0.0541f);
1456  get(dnn::ele_deltaPhiEleClusterTrackAtCalo) =
1457  getValueNorm(electrons.at(index_ele).deltaPhiEleClusterTrackAtCalo(), 0.0002f, 0.0553f);
1458  get(dnn::ele_deltaPhiSuperClusterTrackAtVtx) =
1459  getValueNorm(electrons.at(index_ele).deltaPhiSuperClusterTrackAtVtx(), 0.0001f, 0.0523f);
1460  get(dnn::ele_deltaPhiSeedClusterTrackAtCalo) =
1461  getValueNorm(electrons.at(index_ele).deltaPhiSeedClusterTrackAtCalo(), 0.0004f, 0.0777f);
1462  get(dnn::ele_mvaInput_earlyBrem) = getValue(electrons.at(index_ele).mvaInput().earlyBrem);
1463  get(dnn::ele_mvaInput_lateBrem) = getValue(electrons.at(index_ele).mvaInput().lateBrem);
1464  get(dnn::ele_mvaInput_sigmaEtaEta) =
1465  getValueNorm(electrons.at(index_ele).mvaInput().sigmaEtaEta, 0.0008f, 0.0052f);
1466  get(dnn::ele_mvaInput_hadEnergy) = getValueNorm(electrons.at(index_ele).mvaInput().hadEnergy, 14.04f, 69.48f);
1467  get(dnn::ele_mvaInput_deltaEta) = getValueNorm(electrons.at(index_ele).mvaInput().deltaEta, 0.0099f, 0.0851f);
1468 
1469  const auto& gsfTrack = electrons.at(index_ele).gsfTrack();
1470  if (gsfTrack.isNonnull()) {
1471  get(dnn::ele_gsfTrack_normalizedChi2) = getValueNorm(gsfTrack->normalizedChi2(), 3.049f, 10.39f);
1472  get(dnn::ele_gsfTrack_numberOfValidHits) = getValueNorm(gsfTrack->numberOfValidHits(), 16.52f, 2.806f);
1473  get(dnn::ele_rel_gsfTrack_pt) =
1474  getValueNorm(gsfTrack->pt() / electrons.at(index_ele).polarP4().pt(), 1.355f, 16.81f);
1475  get(dnn::ele_gsfTrack_pt_sig) = getValueNorm(gsfTrack->pt() / gsfTrack->ptError(), 5.046f, 3.119f);
1476  }
1477  const auto& closestCtfTrack = electrons.at(index_ele).closestCtfTrackRef();
1478  const bool has_closestCtfTrack = closestCtfTrack.isNonnull();
1479  if (has_closestCtfTrack) {
1480  get(dnn::ele_has_closestCtfTrack) = has_closestCtfTrack;
1481  get(dnn::ele_closestCtfTrack_normalizedChi2) = getValueNorm(closestCtfTrack->normalizedChi2(), 2.411f, 6.98f);
1482  get(dnn::ele_closestCtfTrack_numberOfValidHits) =
1483  getValueNorm(closestCtfTrack->numberOfValidHits(), 15.16f, 5.26f);
1484  }
1485  }
1486  checkInputs(inputs, is_inner ? "egamma_inner_block" : "egamma_outer_block", dnn::NumberOfInputs);
1487  }
1488 
1490  const TauType& tau,
1491  const reco::Vertex& pv,
1492  double rho,
1493  const pat::MuonCollection& muons,
1494  const pat::PackedCandidateCollection& pfCands,
1495  const Cell& cell_map,
1496  bool is_inner) {
1497  namespace dnn = dnn_inputs_2017_v2::MuonBlockInputs;
1498 
1499  tensorflow::Tensor& inputs = *muonTensor_.at(is_inner);
1500 
1501  const auto& get = [&](int var_index) -> float& { return inputs.tensor<float, 4>()(idx, 0, 0, var_index); };
1502 
1503  const bool valid_index_pf_muon = cell_map.count(CellObjectType::PfCand_muon);
1504  const bool valid_index_muon = cell_map.count(CellObjectType::Muon);
1505 
1506  if (!cell_map.empty()) {
1507  get(dnn::rho) = getValueNorm(rho, 21.49f, 9.713f);
1508  get(dnn::tau_pt) = getValueLinear(tau.polarP4().pt(), 20.f, 1000.f, true);
1509  get(dnn::tau_eta) = getValueLinear(tau.polarP4().eta(), -2.3f, 2.3f, false);
1510  get(dnn::tau_inside_ecal_crack) = getValue(isInEcalCrack(tau.polarP4().eta()));
1511  }
1512  if (valid_index_pf_muon) {
1513  size_t index_pf_muon = cell_map.at(CellObjectType::PfCand_muon);
1514 
1515  get(dnn::pfCand_muon_valid) = valid_index_pf_muon;
1516  get(dnn::pfCand_muon_rel_pt) = getValueNorm(pfCands.at(index_pf_muon).polarP4().pt() / tau.polarP4().pt(),
1517  is_inner ? 0.9509f : 0.0861f,
1518  is_inner ? 0.4294f : 0.4065f);
1519  get(dnn::pfCand_muon_deta) = getValueLinear(pfCands.at(index_pf_muon).polarP4().eta() - tau.polarP4().eta(),
1520  is_inner ? -0.1f : -0.5f,
1521  is_inner ? 0.1f : 0.5f,
1522  false);
1523  get(dnn::pfCand_muon_dphi) = getValueLinear(dPhi(tau.polarP4(), pfCands.at(index_pf_muon).polarP4()),
1524  is_inner ? -0.1f : -0.5f,
1525  is_inner ? 0.1f : 0.5f,
1526  false);
1527  get(dnn::pfCand_muon_pvAssociationQuality) =
1528  getValueLinear<int>(pfCands.at(index_pf_muon).pvAssociationQuality(), 0, 7, true);
1529  get(dnn::pfCand_muon_fromPV) = getValueLinear<int>(pfCands.at(index_pf_muon).fromPV(), 0, 3, true);
1530  get(dnn::pfCand_muon_puppiWeight) = getValue(pfCands.at(index_pf_muon).puppiWeight());
1531  get(dnn::pfCand_muon_charge) = getValue(pfCands.at(index_pf_muon).charge());
1532  get(dnn::pfCand_muon_lostInnerHits) = getValue<int>(pfCands.at(index_pf_muon).lostInnerHits());
1533  get(dnn::pfCand_muon_numberOfPixelHits) =
1534  getValueLinear(pfCands.at(index_pf_muon).numberOfPixelHits(), 0, 11, true);
1535  get(dnn::pfCand_muon_vertex_dx) =
1536  getValueNorm(pfCands.at(index_pf_muon).vertex().x() - pv.position().x(), -0.0007f, 0.6869f);
1537  get(dnn::pfCand_muon_vertex_dy) =
1538  getValueNorm(pfCands.at(index_pf_muon).vertex().y() - pv.position().y(), 0.0001f, 0.6784f);
1539  get(dnn::pfCand_muon_vertex_dz) =
1540  getValueNorm(pfCands.at(index_pf_muon).vertex().z() - pv.position().z(), -0.0117f, 4.097f);
1541  get(dnn::pfCand_muon_vertex_dx_tauFL) = getValueNorm(
1542  pfCands.at(index_pf_muon).vertex().x() - pv.position().x() - tau.flightLength().x(), -0.0001f, 0.8642f);
1543  get(dnn::pfCand_muon_vertex_dy_tauFL) = getValueNorm(
1544  pfCands.at(index_pf_muon).vertex().y() - pv.position().y() - tau.flightLength().y(), 0.0004f, 0.8561f);
1545  get(dnn::pfCand_muon_vertex_dz_tauFL) = getValueNorm(
1546  pfCands.at(index_pf_muon).vertex().z() - pv.position().z() - tau.flightLength().z(), -0.0118f, 4.405f);
1547 
1548  const bool hasTrackDetails = pfCands.at(index_pf_muon).hasTrackDetails();
1549  if (hasTrackDetails) {
1550  get(dnn::pfCand_muon_hasTrackDetails) = hasTrackDetails;
1551  get(dnn::pfCand_muon_dxy) = getValueNorm(pfCands.at(index_pf_muon).dxy(), -0.0045f, 0.9655f);
1552  get(dnn::pfCand_muon_dxy_sig) = getValueNorm(
1553  std::abs(pfCands.at(index_pf_muon).dxy()) / pfCands.at(index_pf_muon).dxyError(), 4.575f, 42.36f);
1554  get(dnn::pfCand_muon_dz) = getValueNorm(pfCands.at(index_pf_muon).dz(), -0.0117f, 4.097f);
1555  get(dnn::pfCand_muon_dz_sig) = getValueNorm(
1556  std::abs(pfCands.at(index_pf_muon).dz()) / pfCands.at(index_pf_muon).dzError(), 80.37f, 343.3f);
1557  get(dnn::pfCand_muon_track_chi2_ndof) = getValueNorm(
1558  pfCands.at(index_pf_muon).pseudoTrack().chi2() / pfCands.at(index_pf_muon).pseudoTrack().ndof(),
1559  0.69f,
1560  1.711f);
1561  get(dnn::pfCand_muon_track_ndof) = getValueNorm(pfCands.at(index_pf_muon).pseudoTrack().ndof(), 17.5f, 5.11f);
1562  }
1563  }
1564  if (valid_index_muon) {
1565  size_t index_muon = cell_map.at(CellObjectType::Muon);
1566 
1567  get(dnn::muon_valid) = valid_index_muon;
1568  get(dnn::muon_rel_pt) = getValueNorm(muons.at(index_muon).polarP4().pt() / tau.polarP4().pt(),
1569  is_inner ? 0.7966f : 0.2678f,
1570  is_inner ? 3.402f : 3.592f);
1571  get(dnn::muon_deta) = getValueLinear(muons.at(index_muon).polarP4().eta() - tau.polarP4().eta(),
1572  is_inner ? -0.1f : -0.5f,
1573  is_inner ? 0.1f : 0.5f,
1574  false);
1575  get(dnn::muon_dphi) = getValueLinear(
1576  dPhi(tau.polarP4(), muons.at(index_muon).polarP4()), is_inner ? -0.1f : -0.5f, is_inner ? 0.1f : 0.5f, false);
1577  get(dnn::muon_dxy) = getValueNorm(muons.at(index_muon).dB(pat::Muon::PV2D), 0.0019f, 1.039f);
1578  get(dnn::muon_dxy_sig) =
1579  getValueNorm(std::abs(muons.at(index_muon).dB(pat::Muon::PV2D)) / muons.at(index_muon).edB(pat::Muon::PV2D),
1580  8.98f,
1581  71.17f);
1582 
1583  const bool normalizedChi2_valid =
1584  muons.at(index_muon).globalTrack().isNonnull() && muons.at(index_muon).normChi2() >= 0;
1585  if (normalizedChi2_valid) {
1586  get(dnn::muon_normalizedChi2_valid) = normalizedChi2_valid;
1587  get(dnn::muon_normalizedChi2) = getValueNorm(muons.at(index_muon).normChi2(), 21.52f, 265.8f);
1588  if (muons.at(index_muon).innerTrack().isNonnull())
1589  get(dnn::muon_numberOfValidHits) = getValueNorm(muons.at(index_muon).numberOfValidHits(), 21.84f, 10.59f);
1590  }
1591  get(dnn::muon_segmentCompatibility) = getValue(muons.at(index_muon).segmentCompatibility());
1592  get(dnn::muon_caloCompatibility) = getValue(muons.at(index_muon).caloCompatibility());
1593 
1594  const bool pfEcalEnergy_valid = muons.at(index_muon).pfEcalEnergy() >= 0;
1595  if (pfEcalEnergy_valid) {
1596  get(dnn::muon_pfEcalEnergy_valid) = pfEcalEnergy_valid;
1597  get(dnn::muon_rel_pfEcalEnergy) =
1598  getValueNorm(muons.at(index_muon).pfEcalEnergy() / muons.at(index_muon).polarP4().pt(), 0.2273f, 0.4865f);
1599  }
1600 
1601  MuonHitMatchV2 hit_match(muons.at(index_muon));
1602  static const std::map<int, std::pair<int, int>> muonMatchHitVars = {
1603  {MuonSubdetId::DT, {dnn::muon_n_matches_DT_1, dnn::muon_n_hits_DT_1}},
1604  {MuonSubdetId::CSC, {dnn::muon_n_matches_CSC_1, dnn::muon_n_hits_CSC_1}},
1605  {MuonSubdetId::RPC, {dnn::muon_n_matches_RPC_1, dnn::muon_n_hits_RPC_1}}};
1606 
1607  static const std::map<int, std::vector<float>> muonMatchVarLimits = {
1608  {MuonSubdetId::DT, {2, 2, 2, 2}}, {MuonSubdetId::CSC, {6, 2, 2, 2}}, {MuonSubdetId::RPC, {7, 6, 4, 4}}};
1609 
1610  static const std::map<int, std::vector<float>> muonHitVarLimits = {{MuonSubdetId::DT, {12, 12, 12, 8}},
1611  {MuonSubdetId::CSC, {24, 12, 12, 12}},
1612  {MuonSubdetId::RPC, {4, 4, 2, 2}}};
1613 
1614  for (int subdet : hit_match.MuonHitMatchV2::consideredSubdets()) {
1615  const auto& matchHitVar = muonMatchHitVars.at(subdet);
1616  const auto& matchLimits = muonMatchVarLimits.at(subdet);
1617  const auto& hitLimits = muonHitVarLimits.at(subdet);
1618  for (int station = MuonHitMatchV2::first_station_id; station <= MuonHitMatchV2::last_station_id; ++station) {
1619  const unsigned n_matches = hit_match.nMatches(subdet, station);
1620  const unsigned n_hits = hit_match.nHits(subdet, station);
1621  get(matchHitVar.first + station - 1) = getValueLinear(n_matches, 0, matchLimits.at(station - 1), true);
1622  get(matchHitVar.second + station - 1) = getValueLinear(n_hits, 0, hitLimits.at(station - 1), true);
1623  }
1624  }
1625  }
1626  checkInputs(inputs, is_inner ? "muon_inner_block" : "muon_outer_block", dnn::NumberOfInputs);
1627  }
1628 
1630  const TauType& tau,
1631  const reco::Vertex& pv,
1632  double rho,
1633  const pat::PackedCandidateCollection& pfCands,
1634  const Cell& cell_map,
1635  bool is_inner) {
1636  namespace dnn = dnn_inputs_2017_v2::HadronBlockInputs;
1637 
1638  tensorflow::Tensor& inputs = *hadronsTensor_.at(is_inner);
1639 
1640  const auto& get = [&](int var_index) -> float& { return inputs.tensor<float, 4>()(idx, 0, 0, var_index); };
1641 
1642  const bool valid_chH = cell_map.count(CellObjectType::PfCand_chargedHadron);
1643  const bool valid_nH = cell_map.count(CellObjectType::PfCand_neutralHadron);
1644 
1645  if (!cell_map.empty()) {
1646  get(dnn::rho) = getValueNorm(rho, 21.49f, 9.713f);
1647  get(dnn::tau_pt) = getValueLinear(tau.polarP4().pt(), 20.f, 1000.f, true);
1648  get(dnn::tau_eta) = getValueLinear(tau.polarP4().eta(), -2.3f, 2.3f, false);
1649  get(dnn::tau_inside_ecal_crack) = getValue(isInEcalCrack(tau.polarP4().eta()));
1650  }
1651  if (valid_chH) {
1652  size_t index_chH = cell_map.at(CellObjectType::PfCand_chargedHadron);
1653 
1654  get(dnn::pfCand_chHad_valid) = valid_chH;
1655  get(dnn::pfCand_chHad_rel_pt) = getValueNorm(pfCands.at(index_chH).polarP4().pt() / tau.polarP4().pt(),
1656  is_inner ? 0.2564f : 0.0194f,
1657  is_inner ? 0.8607f : 0.1865f);
1658  get(dnn::pfCand_chHad_deta) = getValueLinear(pfCands.at(index_chH).polarP4().eta() - tau.polarP4().eta(),
1659  is_inner ? -0.1f : -0.5f,
1660  is_inner ? 0.1f : 0.5f,
1661  false);
1662  get(dnn::pfCand_chHad_dphi) = getValueLinear(dPhi(tau.polarP4(), pfCands.at(index_chH).polarP4()),
1663  is_inner ? -0.1f : -0.5f,
1664  is_inner ? 0.1f : 0.5f,
1665  false);
1666  get(dnn::pfCand_chHad_leadChargedHadrCand) = getValue(
1667  &pfCands.at(index_chH) == dynamic_cast<const pat::PackedCandidate*>(tau.leadChargedHadrCand().get()));
1668  get(dnn::pfCand_chHad_pvAssociationQuality) =
1669  getValueLinear<int>(pfCands.at(index_chH).pvAssociationQuality(), 0, 7, true);
1670  get(dnn::pfCand_chHad_fromPV) = getValueLinear<int>(pfCands.at(index_chH).fromPV(), 0, 3, true);
1671  get(dnn::pfCand_chHad_puppiWeight) = getValue(pfCands.at(index_chH).puppiWeight());
1672  get(dnn::pfCand_chHad_puppiWeightNoLep) = getValue(pfCands.at(index_chH).puppiWeightNoLep());
1673  get(dnn::pfCand_chHad_charge) = getValue(pfCands.at(index_chH).charge());
1674  get(dnn::pfCand_chHad_lostInnerHits) = getValue<int>(pfCands.at(index_chH).lostInnerHits());
1675  get(dnn::pfCand_chHad_numberOfPixelHits) = getValueLinear(pfCands.at(index_chH).numberOfPixelHits(), 0, 12, true);
1676  get(dnn::pfCand_chHad_vertex_dx) =
1677  getValueNorm(pfCands.at(index_chH).vertex().x() - pv.position().x(), 0.0005f, 1.735f);
1678  get(dnn::pfCand_chHad_vertex_dy) =
1679  getValueNorm(pfCands.at(index_chH).vertex().y() - pv.position().y(), -0.0008f, 1.752f);
1680  get(dnn::pfCand_chHad_vertex_dz) =
1681  getValueNorm(pfCands.at(index_chH).vertex().z() - pv.position().z(), -0.0201f, 8.333f);
1682  get(dnn::pfCand_chHad_vertex_dx_tauFL) = getValueNorm(
1683  pfCands.at(index_chH).vertex().x() - pv.position().x() - tau.flightLength().x(), -0.0014f, 1.93f);
1684  get(dnn::pfCand_chHad_vertex_dy_tauFL) = getValueNorm(
1685  pfCands.at(index_chH).vertex().y() - pv.position().y() - tau.flightLength().y(), 0.0022f, 1.948f);
1686  get(dnn::pfCand_chHad_vertex_dz_tauFL) = getValueNorm(
1687  pfCands.at(index_chH).vertex().z() - pv.position().z() - tau.flightLength().z(), -0.0138f, 8.622f);
1688 
1689  const bool hasTrackDetails = pfCands.at(index_chH).hasTrackDetails();
1690  if (hasTrackDetails) {
1691  get(dnn::pfCand_chHad_hasTrackDetails) = hasTrackDetails;
1692  get(dnn::pfCand_chHad_dxy) = getValueNorm(pfCands.at(index_chH).dxy(), -0.012f, 2.386f);
1693  get(dnn::pfCand_chHad_dxy_sig) =
1694  getValueNorm(std::abs(pfCands.at(index_chH).dxy()) / pfCands.at(index_chH).dxyError(), 6.417f, 36.28f);
1695  get(dnn::pfCand_chHad_dz) = getValueNorm(pfCands.at(index_chH).dz(), -0.0246f, 7.618f);
1696  get(dnn::pfCand_chHad_dz_sig) =
1697  getValueNorm(std::abs(pfCands.at(index_chH).dz()) / pfCands.at(index_chH).dzError(), 301.3f, 491.1f);
1698  get(dnn::pfCand_chHad_track_chi2_ndof) =
1699  pfCands.at(index_chH).pseudoTrack().ndof() > 0
1700  ? getValueNorm(pfCands.at(index_chH).pseudoTrack().chi2() / pfCands.at(index_chH).pseudoTrack().ndof(),
1701  0.7876f,
1702  3.694f)
1703  : 0;
1704  get(dnn::pfCand_chHad_track_ndof) =
1705  pfCands.at(index_chH).pseudoTrack().ndof() > 0
1706  ? getValueNorm(pfCands.at(index_chH).pseudoTrack().ndof(), 13.92f, 6.581f)
1707  : 0;
1708  }
1709  float hcal_fraction = 0.;
1710  if (pfCands.at(index_chH).pdgId() == 1 || pfCands.at(index_chH).pdgId() == 130) {
1711  hcal_fraction = pfCands.at(index_chH).hcalFraction();
1712  } else if (pfCands.at(index_chH).isIsolatedChargedHadron()) {
1713  hcal_fraction = pfCands.at(index_chH).rawHcalFraction();
1714  }
1715  get(dnn::pfCand_chHad_hcalFraction) = getValue(hcal_fraction);
1716  get(dnn::pfCand_chHad_rawCaloFraction) = getValueLinear(pfCands.at(index_chH).rawCaloFraction(), 0.f, 2.6f, true);
1717  }
1718  if (valid_nH) {
1719  size_t index_nH = cell_map.at(CellObjectType::PfCand_neutralHadron);
1720 
1721  get(dnn::pfCand_nHad_valid) = valid_nH;
1722  get(dnn::pfCand_nHad_rel_pt) = getValueNorm(pfCands.at(index_nH).polarP4().pt() / tau.polarP4().pt(),
1723  is_inner ? 0.3163f : 0.0502f,
1724  is_inner ? 0.2769f : 0.4266f);
1725  get(dnn::pfCand_nHad_deta) = getValueLinear(pfCands.at(index_nH).polarP4().eta() - tau.polarP4().eta(),
1726  is_inner ? -0.1f : -0.5f,
1727  is_inner ? 0.1f : 0.5f,
1728  false);
1729  get(dnn::pfCand_nHad_dphi) = getValueLinear(
1730  dPhi(tau.polarP4(), pfCands.at(index_nH).polarP4()), is_inner ? -0.1f : -0.5f, is_inner ? 0.1f : 0.5f, false);
1731  get(dnn::pfCand_nHad_puppiWeight) = getValue(pfCands.at(index_nH).puppiWeight());
1732  get(dnn::pfCand_nHad_puppiWeightNoLep) = getValue(pfCands.at(index_nH).puppiWeightNoLep());
1733  float hcal_fraction = 0.;
1734  if (pfCands.at(index_nH).pdgId() == 1 || pfCands.at(index_nH).pdgId() == 130) {
1735  hcal_fraction = pfCands.at(index_nH).hcalFraction();
1736  } else if (pfCands.at(index_nH).isIsolatedChargedHadron()) {
1737  hcal_fraction = pfCands.at(index_nH).rawHcalFraction();
1738  }
1739  get(dnn::pfCand_nHad_hcalFraction) = getValue(hcal_fraction);
1740  }
1741  checkInputs(inputs, is_inner ? "hadron_inner_block" : "hadron_outer_block", dnn::NumberOfInputs);
1742  }
1743 
1744  template <typename dnn>
1745  tensorflow::Tensor createInputsV1(const TauType& tau,
1747  const MuonCollection& muons) const {
1748  static constexpr bool check_all_set = false;
1749  static constexpr float default_value_for_set_check = -42;
1750 
1751  tensorflow::Tensor inputs(tensorflow::DT_FLOAT, {1, dnn_inputs_2017v1::NumberOfInputs});
1752  const auto& get = [&](int var_index) -> float& { return inputs.matrix<float>()(0, var_index); };
1753  auto leadChargedHadrCand = dynamic_cast<const pat::PackedCandidate*>(tau.leadChargedHadrCand().get());
1754 
1755  if (check_all_set) {
1756  for (int var_index = 0; var_index < dnn::NumberOfInputs; ++var_index) {
1757  get(var_index) = default_value_for_set_check;
1758  }
1759  }
1760 
1761  get(dnn::pt) = tau.p4().pt();
1762  get(dnn::eta) = tau.p4().eta();
1763  get(dnn::mass) = tau.p4().mass();
1764  get(dnn::decayMode) = tau.decayMode();
1765  get(dnn::chargedIsoPtSum) = tau.tauID("chargedIsoPtSum");
1766  get(dnn::neutralIsoPtSum) = tau.tauID("neutralIsoPtSum");
1767  get(dnn::neutralIsoPtSumWeight) = tau.tauID("neutralIsoPtSumWeight");
1768  get(dnn::photonPtSumOutsideSignalCone) = tau.tauID("photonPtSumOutsideSignalCone");
1769  get(dnn::puCorrPtSum) = tau.tauID("puCorrPtSum");
1770  get(dnn::dxy) = tau.dxy();
1771  get(dnn::dxy_sig) = tau.dxy_Sig();
1772  get(dnn::dz) = leadChargedHadrCand ? leadChargedHadrCand->dz() : default_value;
1773  get(dnn::ip3d) = tau.ip3d();
1774  get(dnn::ip3d_sig) = tau.ip3d_Sig();
1775  get(dnn::hasSecondaryVertex) = tau.hasSecondaryVertex();
1776  get(dnn::flightLength_r) = tau.flightLength().R();
1777  get(dnn::flightLength_dEta) = dEta(tau.flightLength(), tau.p4());
1778  get(dnn::flightLength_dPhi) = dPhi(tau.flightLength(), tau.p4());
1779  get(dnn::flightLength_sig) = tau.flightLengthSig();
1780  get(dnn::leadChargedHadrCand_pt) = leadChargedHadrCand ? leadChargedHadrCand->p4().Pt() : default_value;
1781  get(dnn::leadChargedHadrCand_dEta) =
1782  leadChargedHadrCand ? dEta(leadChargedHadrCand->p4(), tau.p4()) : default_value;
1783  get(dnn::leadChargedHadrCand_dPhi) =
1784  leadChargedHadrCand ? dPhi(leadChargedHadrCand->p4(), tau.p4()) : default_value;
1785  get(dnn::leadChargedHadrCand_mass) = leadChargedHadrCand ? leadChargedHadrCand->p4().mass() : default_value;
1790  get(dnn::leadingTrackNormChi2) = tau.leadingTrackNormChi2();
1791  get(dnn::e_ratio) = reco::tau::eratio(tau);
1792  get(dnn::gj_angle_diff) = calculateGottfriedJacksonAngleDifference(tau);
1793  get(dnn::n_photons) = reco::tau::n_photons_total(tau);
1794  get(dnn::emFraction) = tau.emFraction_MVA();
1795  get(dnn::has_gsf_track) = leadChargedHadrCand && std::abs(leadChargedHadrCand->pdgId()) == 11;
1796  get(dnn::inside_ecal_crack) = isInEcalCrack(tau.p4().Eta());
1797  auto gsf_ele = findMatchedElectron(tau, electrons, 0.3);
1798  get(dnn::gsf_ele_matched) = gsf_ele != nullptr;
1799  get(dnn::gsf_ele_pt) = gsf_ele != nullptr ? gsf_ele->p4().Pt() : default_value;
1800  get(dnn::gsf_ele_dEta) = gsf_ele != nullptr ? dEta(gsf_ele->p4(), tau.p4()) : default_value;
1801  get(dnn::gsf_ele_dPhi) = gsf_ele != nullptr ? dPhi(gsf_ele->p4(), tau.p4()) : default_value;
1802  get(dnn::gsf_ele_mass) = gsf_ele != nullptr ? gsf_ele->p4().mass() : default_value;
1803  calculateElectronClusterVars(gsf_ele, get(dnn::gsf_ele_Ee), get(dnn::gsf_ele_Egamma));
1804  get(dnn::gsf_ele_Pin) = gsf_ele != nullptr ? gsf_ele->trackMomentumAtVtx().R() : default_value;
1805  get(dnn::gsf_ele_Pout) = gsf_ele != nullptr ? gsf_ele->trackMomentumOut().R() : default_value;
1806  get(dnn::gsf_ele_EtotOverPin) = get(dnn::gsf_ele_Pin) > 0
1807  ? (get(dnn::gsf_ele_Ee) + get(dnn::gsf_ele_Egamma)) / get(dnn::gsf_ele_Pin)
1808  : default_value;
1809  get(dnn::gsf_ele_Eecal) = gsf_ele != nullptr ? gsf_ele->ecalEnergy() : default_value;
1810  get(dnn::gsf_ele_dEta_SeedClusterTrackAtCalo) =
1811  gsf_ele != nullptr ? gsf_ele->deltaEtaSeedClusterTrackAtCalo() : default_value;
1812  get(dnn::gsf_ele_dPhi_SeedClusterTrackAtCalo) =
1813  gsf_ele != nullptr ? gsf_ele->deltaPhiSeedClusterTrackAtCalo() : default_value;
1814  get(dnn::gsf_ele_mvaIn_sigmaEtaEta) = gsf_ele != nullptr ? gsf_ele->mvaInput().sigmaEtaEta : default_value;
1815  get(dnn::gsf_ele_mvaIn_hadEnergy) = gsf_ele != nullptr ? gsf_ele->mvaInput().hadEnergy : default_value;
1816  get(dnn::gsf_ele_mvaIn_deltaEta) = gsf_ele != nullptr ? gsf_ele->mvaInput().deltaEta : default_value;
1817 
1818  get(dnn::gsf_ele_Chi2NormGSF) = default_value;
1819  get(dnn::gsf_ele_GSFNumHits) = default_value;
1820  get(dnn::gsf_ele_GSFTrackResol) = default_value;
1821  get(dnn::gsf_ele_GSFTracklnPt) = default_value;
1822  if (gsf_ele != nullptr && gsf_ele->gsfTrack().isNonnull()) {
1823  get(dnn::gsf_ele_Chi2NormGSF) = gsf_ele->gsfTrack()->normalizedChi2();
1824  get(dnn::gsf_ele_GSFNumHits) = gsf_ele->gsfTrack()->numberOfValidHits();
1825  if (gsf_ele->gsfTrack()->pt() > 0) {
1826  get(dnn::gsf_ele_GSFTrackResol) = gsf_ele->gsfTrack()->ptError() / gsf_ele->gsfTrack()->pt();
1827  get(dnn::gsf_ele_GSFTracklnPt) = std::log10(gsf_ele->gsfTrack()->pt());
1828  }
1829  }
1830 
1831  get(dnn::gsf_ele_Chi2NormKF) = default_value;
1832  get(dnn::gsf_ele_KFNumHits) = default_value;
1833  if (gsf_ele != nullptr && gsf_ele->closestCtfTrackRef().isNonnull()) {
1834  get(dnn::gsf_ele_Chi2NormKF) = gsf_ele->closestCtfTrackRef()->normalizedChi2();
1835  get(dnn::gsf_ele_KFNumHits) = gsf_ele->closestCtfTrackRef()->numberOfValidHits();
1836  }
1837  get(dnn::leadChargedCand_etaAtEcalEntrance) = tau.etaAtEcalEntranceLeadChargedCand();
1838  get(dnn::leadChargedCand_pt) = tau.ptLeadChargedCand();
1839 
1840  get(dnn::leadChargedHadrCand_HoP) = default_value;
1841  get(dnn::leadChargedHadrCand_EoP) = default_value;
1842  if (tau.leadChargedHadrCand()->pt() > 0) {
1843  get(dnn::leadChargedHadrCand_HoP) = tau.hcalEnergyLeadChargedHadrCand() / tau.leadChargedHadrCand()->pt();
1844  get(dnn::leadChargedHadrCand_EoP) = tau.ecalEnergyLeadChargedHadrCand() / tau.leadChargedHadrCand()->pt();
1845  }
1846 
1847  MuonHitMatchV1 muon_hit_match;
1848  if (tau.leadPFChargedHadrCand().isNonnull() && tau.leadPFChargedHadrCand()->muonRef().isNonnull())
1849  muon_hit_match.addMatchedMuon(*tau.leadPFChargedHadrCand()->muonRef(), tau);
1850 
1851  auto matched_muons = muon_hit_match.findMatchedMuons(tau, muons, 0.3, 5);
1852  for (auto muon : matched_muons)
1853  muon_hit_match.addMatchedMuon(*muon, tau);
1854  muon_hit_match.fillTensor<dnn>(get, tau, default_value);
1855 
1856  LorentzVectorXYZ signalChargedHadrCands_sumIn, signalChargedHadrCands_sumOut;
1858  tau.signalChargedHadrCands(),
1859  signalChargedHadrCands_sumIn,
1860  signalChargedHadrCands_sumOut,
1861  get(dnn::signalChargedHadrCands_sum_innerSigCone_pt),
1862  get(dnn::signalChargedHadrCands_sum_innerSigCone_dEta),
1863  get(dnn::signalChargedHadrCands_sum_innerSigCone_dPhi),
1864  get(dnn::signalChargedHadrCands_sum_innerSigCone_mass),
1865  get(dnn::signalChargedHadrCands_sum_outerSigCone_pt),
1866  get(dnn::signalChargedHadrCands_sum_outerSigCone_dEta),
1867  get(dnn::signalChargedHadrCands_sum_outerSigCone_dPhi),
1868  get(dnn::signalChargedHadrCands_sum_outerSigCone_mass),
1869  get(dnn::signalChargedHadrCands_nTotal_innerSigCone),
1870  get(dnn::signalChargedHadrCands_nTotal_outerSigCone));
1871 
1872  LorentzVectorXYZ signalNeutrHadrCands_sumIn, signalNeutrHadrCands_sumOut;
1874  tau.signalNeutrHadrCands(),
1875  signalNeutrHadrCands_sumIn,
1876  signalNeutrHadrCands_sumOut,
1877  get(dnn::signalNeutrHadrCands_sum_innerSigCone_pt),
1878  get(dnn::signalNeutrHadrCands_sum_innerSigCone_dEta),
1879  get(dnn::signalNeutrHadrCands_sum_innerSigCone_dPhi),
1880  get(dnn::signalNeutrHadrCands_sum_innerSigCone_mass),
1881  get(dnn::signalNeutrHadrCands_sum_outerSigCone_pt),
1882  get(dnn::signalNeutrHadrCands_sum_outerSigCone_dEta),
1883  get(dnn::signalNeutrHadrCands_sum_outerSigCone_dPhi),
1884  get(dnn::signalNeutrHadrCands_sum_outerSigCone_mass),
1885  get(dnn::signalNeutrHadrCands_nTotal_innerSigCone),
1886  get(dnn::signalNeutrHadrCands_nTotal_outerSigCone));
1887 
1888  LorentzVectorXYZ signalGammaCands_sumIn, signalGammaCands_sumOut;
1890  tau.signalGammaCands(),
1891  signalGammaCands_sumIn,
1892  signalGammaCands_sumOut,
1893  get(dnn::signalGammaCands_sum_innerSigCone_pt),
1894  get(dnn::signalGammaCands_sum_innerSigCone_dEta),
1895  get(dnn::signalGammaCands_sum_innerSigCone_dPhi),
1896  get(dnn::signalGammaCands_sum_innerSigCone_mass),
1897  get(dnn::signalGammaCands_sum_outerSigCone_pt),
1898  get(dnn::signalGammaCands_sum_outerSigCone_dEta),
1899  get(dnn::signalGammaCands_sum_outerSigCone_dPhi),
1900  get(dnn::signalGammaCands_sum_outerSigCone_mass),
1901  get(dnn::signalGammaCands_nTotal_innerSigCone),
1902  get(dnn::signalGammaCands_nTotal_outerSigCone));
1903 
1904  LorentzVectorXYZ isolationChargedHadrCands_sum;
1906  tau.isolationChargedHadrCands(),
1907  isolationChargedHadrCands_sum,
1908  get(dnn::isolationChargedHadrCands_sum_pt),
1909  get(dnn::isolationChargedHadrCands_sum_dEta),
1910  get(dnn::isolationChargedHadrCands_sum_dPhi),
1911  get(dnn::isolationChargedHadrCands_sum_mass),
1912  get(dnn::isolationChargedHadrCands_nTotal));
1913 
1914  LorentzVectorXYZ isolationNeutrHadrCands_sum;
1916  tau.isolationNeutrHadrCands(),
1917  isolationNeutrHadrCands_sum,
1918  get(dnn::isolationNeutrHadrCands_sum_pt),
1919  get(dnn::isolationNeutrHadrCands_sum_dEta),
1920  get(dnn::isolationNeutrHadrCands_sum_dPhi),
1921  get(dnn::isolationNeutrHadrCands_sum_mass),
1922  get(dnn::isolationNeutrHadrCands_nTotal));
1923 
1924  LorentzVectorXYZ isolationGammaCands_sum;
1926  tau.isolationGammaCands(),
1927  isolationGammaCands_sum,
1928  get(dnn::isolationGammaCands_sum_pt),
1929  get(dnn::isolationGammaCands_sum_dEta),
1930  get(dnn::isolationGammaCands_sum_dPhi),
1931  get(dnn::isolationGammaCands_sum_mass),
1932  get(dnn::isolationGammaCands_nTotal));
1933 
1934  get(dnn::tau_visMass_innerSigCone) = (signalGammaCands_sumIn + signalChargedHadrCands_sumIn).mass();
1935 
1936  if (check_all_set) {
1937  for (int var_index = 0; var_index < dnn::NumberOfInputs; ++var_index) {
1938  if (get(var_index) == default_value_for_set_check)
1939  throw cms::Exception("DeepTauId: variable with index = ") << var_index << " is not set.";
1940  }
1941  }
1942 
1943  return inputs;
1944  }
1945 
1946  static void calculateElectronClusterVars(const pat::Electron* ele, float& elecEe, float& elecEgamma) {
1947  if (ele) {
1948  elecEe = elecEgamma = 0;
1949  auto superCluster = ele->superCluster();
1950  if (superCluster.isNonnull() && superCluster.isAvailable() && superCluster->clusters().isNonnull() &&
1951  superCluster->clusters().isAvailable()) {
1952  for (auto iter = superCluster->clustersBegin(); iter != superCluster->clustersEnd(); ++iter) {
1953  const double energy = (*iter)->energy();
1954  if (iter == superCluster->clustersBegin())
1955  elecEe += energy;
1956  else
1957  elecEgamma += energy;
1958  }
1959  }
1960  } else {
1961  elecEe = elecEgamma = default_value;
1962  }
1963  }
1964 
1965  template <typename CandidateCollection>
1968  LorentzVectorXYZ& p4_inner,
1969  LorentzVectorXYZ& p4_outer,
1970  float& pt_inner,
1971  float& dEta_inner,
1972  float& dPhi_inner,
1973  float& m_inner,
1974  float& pt_outer,
1975  float& dEta_outer,
1976  float& dPhi_outer,
1977  float& m_outer,
1978  float& n_inner,
1979  float& n_outer) {
1980  p4_inner = LorentzVectorXYZ(0, 0, 0, 0);
1981  p4_outer = LorentzVectorXYZ(0, 0, 0, 0);
1982  n_inner = 0;
1983  n_outer = 0;
1984 
1985  const double innerSigCone_radius = getInnerSignalConeRadius(tau.pt());
1986  for (const auto& cand : candidates) {
1987  const double dR = reco::deltaR(cand->p4(), tau.leadChargedHadrCand()->p4());
1988  const bool isInside_innerSigCone = dR < innerSigCone_radius;
1989  if (isInside_innerSigCone) {
1990  p4_inner += cand->p4();
1991  ++n_inner;
1992  } else {
1993  p4_outer += cand->p4();
1994  ++n_outer;
1995  }
1996  }
1997 
1998  pt_inner = n_inner != 0 ? p4_inner.Pt() : default_value;
1999  dEta_inner = n_inner != 0 ? dEta(p4_inner, tau.p4()) : default_value;
2000  dPhi_inner = n_inner != 0 ? dPhi(p4_inner, tau.p4()) : default_value;
2001  m_inner = n_inner != 0 ? p4_inner.mass() : default_value;
2002 
2003  pt_outer = n_outer != 0 ? p4_outer.Pt() : default_value;
2004  dEta_outer = n_outer != 0 ? dEta(p4_outer, tau.p4()) : default_value;
2005  dPhi_outer = n_outer != 0 ? dPhi(p4_outer, tau.p4()) : default_value;
2006  m_outer = n_outer != 0 ? p4_outer.mass() : default_value;
2007  }
2008 
2009  template <typename CandidateCollection>
2013  float& pt,
2014  float& d_eta,
2015  float& d_phi,
2016  float& m,
2017  float& n) {
2018  p4 = LorentzVectorXYZ(0, 0, 0, 0);
2019  n = 0;
2020 
2021  for (const auto& cand : candidates) {
2022  p4 += cand->p4();
2023  ++n;
2024  }
2025 
2026  pt = n != 0 ? p4.Pt() : default_value;
2027  d_eta = n != 0 ? dEta(p4, tau.p4()) : default_value;
2028  d_phi = n != 0 ? dPhi(p4, tau.p4()) : default_value;
2029  m = n != 0 ? p4.mass() : default_value;
2030  }
2031 
2032  static double getInnerSignalConeRadius(double pt) {
2033  static constexpr double min_pt = 30., min_radius = 0.05, cone_opening_coef = 3.;
2034  // This is equivalent of the original formula (std::max(std::min(0.1, 3.0/pt), 0.05)
2035  return std::max(cone_opening_coef / std::max(pt, min_pt), min_radius);
2036  }
2037 
2038  // Copied from https://github.com/cms-sw/cmssw/blob/CMSSW_9_4_X/RecoTauTag/RecoTau/plugins/PATTauDiscriminationByMVAIsolationRun2.cc#L218
2039  static bool calculateGottfriedJacksonAngleDifference(const pat::Tau& tau, double& gj_diff) {
2040  if (tau.hasSecondaryVertex()) {
2041  static constexpr double mTau = 1.77682;
2042  const double mAOne = tau.p4().M();
2043  const double pAOneMag = tau.p();
2044  const double argumentThetaGJmax = (std::pow(mTau, 2) - std::pow(mAOne, 2)) / (2 * mTau * pAOneMag);
2045  const double argumentThetaGJmeasured =
2046  tau.p4().Vect().Dot(tau.flightLength()) / (pAOneMag * tau.flightLength().R());
2047  if (std::abs(argumentThetaGJmax) <= 1. && std::abs(argumentThetaGJmeasured) <= 1.) {
2048  double thetaGJmax = std::asin(argumentThetaGJmax);
2049  double thetaGJmeasured = std::acos(argumentThetaGJmeasured);
2050  gj_diff = thetaGJmeasured - thetaGJmax;
2051  return true;
2052  }
2053  }
2054  return false;
2055  }
2056 
2058  double gj_diff;
2060  return static_cast<float>(gj_diff);
2061  return default_value;
2062  }
2063 
2064  static bool isInEcalCrack(double eta) {
2065  const double abs_eta = std::abs(eta);
2066  return abs_eta > 1.46 && abs_eta < 1.558;
2067  }
2068 
2071  double deltaR) {
2072  const double dR2 = deltaR * deltaR;
2073  const pat::Electron* matched_ele = nullptr;
2074  for (const auto& ele : electrons) {
2075  if (reco::deltaR2(tau.p4(), ele.p4()) < dR2 && (!matched_ele || matched_ele->pt() < ele.pt())) {
2076  matched_ele = &ele;
2077  }
2078  }
2079  return matched_ele;
2080  }
2081 
2082 private:
2087  const unsigned version;
2088  const int debug_level;
2089  const bool disable_dxy_pca_;
2090  std::unique_ptr<tensorflow::Tensor> tauBlockTensor_;
2091  std::array<std::unique_ptr<tensorflow::Tensor>, 2> eGammaTensor_, muonTensor_, hadronsTensor_, convTensor_,
2093 };
2094 
reco::tau::pt_weighted_deta_strip
float pt_weighted_deta_strip(const reco::PFTau &tau, int dm)
Definition: PFRecoTauClusterVariables.h:36
deep_tau::DeepTauBase::cache_
const DeepTauCache * cache_
Definition: DeepTauBase.h:104
PDWG_BPHSkim_cff.muons
muons
Definition: PDWG_BPHSkim_cff.py:47
alignBH_cfg.fixed
fixed
Definition: alignBH_cfg.py:54
MuonSubdetId::CSC
static constexpr int CSC
Definition: MuonSubdetId.h:12
DeepTauId::rho_token_
edm::EDGetTokenT< double > rho_token_
Definition: DeepTauId.cc:2085
dnn_inputs_2017_v2
Definition: DeepTauId.cc:153
DeepTauId::calculateElectronClusterVars
static void calculateElectronClusterVars(const pat::Electron *ele, float &elecEe, float &elecEgamma)
Definition: DeepTauId.cc:1946
DeepTauBase.h
DeepTauId::DeepTauId
DeepTauId(const edm::ParameterSet &cfg, const deep_tau::DeepTauCache *cache)
Definition: DeepTauId.cc:844
electrons_cff.bool
bool
Definition: electrons_cff.py:372
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
DeepTauId::getInnerSignalConeRadius
static double getInnerSignalConeRadius(double pt)
Definition: DeepTauId.cc:2032
input
static const std::string input
Definition: EdmProvDump.cc:48
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
deep_tau::DeepTauCache::getSession
tensorflow::Session & getSession(const std::string &name="") const
Definition: DeepTauBase.h:49
SiStripPI::mean
Definition: SiStripPayloadInspectorHelper.h:169
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
DeepTauId::processSignalPFComponents
static void processSignalPFComponents(const pat::Tau &tau, const CandidateCollection &candidates, LorentzVectorXYZ &p4_inner, LorentzVectorXYZ &p4_outer, float &pt_inner, float &dEta_inner, float &dPhi_inner, float &m_inner, float &pt_outer, float &dEta_outer, float &dPhi_outer, float &m_outer, float &n_inner, float &n_outer)
Definition: DeepTauId.cc:1966
runTauDisplay.tau_eta
tau_eta
Definition: runTauDisplay.py:126
dnn_inputs_2017_v2::MuonBlockInputs
Definition: DeepTauId.cc:302
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
metsig::tau
Definition: SignAlgoResolutions.h:49
muon
Definition: MuonCocktails.h:17
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
DeepTauId::isInEcalCrack
static bool isInEcalCrack(double eta)
Definition: DeepTauId.cc:2064
pat::ElectronCollection
std::vector< Electron > ElectronCollection
Definition: Electron.h:36
DeepTauId::pi
static constexpr float pi
Definition: DeepTauId.cc:899
DeepTauId::getValue
static float getValue(T value)
Definition: DeepTauId.cc:902
DeepTauId::setCellConvFeatures
void setCellConvFeatures(tensorflow::Tensor &convTensor, const tensorflow::Tensor &features, unsigned batch_idx, int eta_index, int phi_index)
Definition: DeepTauId.cc:1162
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
DeepTauId::createTauBlockInputs
void createTauBlockInputs(const TauType &tau, const reco::Vertex &pv, double rho)
Definition: DeepTauId.cc:1171
reco::deltaPhi
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
postprocess-scan-build.features
features
Definition: postprocess-scan-build.py:8
edm::EDGetTokenT< ElectronCollection >
Tau3MuMonitor_cff.taus
taus
Definition: Tau3MuMonitor_cff.py:7
relativeConstraints.station
station
Definition: relativeConstraints.py:67
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::tau::eratio
float eratio(const reco::PFTau &tau)
return ratio of energy in ECAL over sum of energy in ECAL and HCAL
Definition: PFRecoTauClusterVariables.cc:78
DeepTauId::findMatchedElectron
static const pat::Electron * findMatchedElectron(const pat::Tau &tau, const pat::ElectronCollection &electrons, double deltaR)
Definition: DeepTauId.cc:2069
fireworks::subdets
static const std::string subdets[7]
Definition: TrackUtils.cc:60
deep_tau::DeepTauBase::pfcandToken_
edm::EDGetTokenT< pat::PackedCandidateCollection > pfcandToken_
Definition: DeepTauBase.h:100
pat::Tau
Analysis-level tau class.
Definition: Tau.h:53
deep_tau::DeepTauBase::OutputCollection
std::map< std::string, Output > OutputCollection
Definition: DeepTauBase.h:82
gather_cfg.cout
cout
Definition: gather_cfg.py:144
objects
Definition: __init__.py:1
ZElectronSkim_cff.rho
rho
Definition: ZElectronSkim_cff.py:38
DeepTauId::hadronsTensor_
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > hadronsTensor_
Definition: DeepTauId.cc:2091
DeepTauId::createInputsV1
tensorflow::Tensor createInputsV1(const TauType &tau, const ElectronCollection &electrons, const MuonCollection &muons) const
Definition: DeepTauId.cc:1745
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
DeepTauId::debug_level
const int debug_level
Definition: DeepTauId.cc:2088
DeepTauId::getValueNorm
static float getValueNorm(T value, float mean, float sigma, float n_sigmas_max=5)
Definition: DeepTauId.cc:917
getRunAppsInfo.grid
grid
Definition: getRunAppsInfo.py:92
deep_tau::DeepTauBase::vtxToken_
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
Definition: DeepTauBase.h:101
DeepTauId::initializeGlobalCache
static std::unique_ptr< deep_tau::DeepTauCache > initializeGlobalCache(const edm::ParameterSet &cfg)
Definition: DeepTauId.cc:892
deep_tau
Definition: DeepTauBase.h:28
DDAxes::x
deep_tau::DeepTauBase::LorentzVectorXYZ
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > LorentzVectorXYZ
Definition: DeepTauBase.h:67
reco::LeafCandidate::pt
double pt() const final
transverse momentum
Definition: LeafCandidate.h:146
pat::Muon
Analysis-level muon class.
Definition: Muon.h:51
DeepTauId::disable_dxy_pca_
const bool disable_dxy_pca_
Definition: DeepTauId.cc:2089
HLT_2018_cff.muon
muon
Definition: HLT_2018_cff.py:10349
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
edm::Handle< TauCollection >
training_settings.idx
idx
Definition: training_settings.py:16
DeepTauId::zeroOutputTensor_
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > zeroOutputTensor_
Definition: DeepTauId.cc:2091
deep_tau::NumberOfOutputs
constexpr int NumberOfOutputs
Definition: DeepTauId.cc:13
cosmictrackSelector_cfi.min_pt
min_pt
Definition: cosmictrackSelector_cfi.py:18
HLT_2018_cff.dEta
dEta
Definition: HLT_2018_cff.py:12289
end
#define end
Definition: vmac.h:39
reco::Muon
Definition: Muon.h:27
DeepTauId::default_value
static constexpr float default_value
Definition: DeepTauId.cc:810
deep_tau::DeepTauBase
Definition: DeepTauBase.h:58
deep_tau::DeepTauBase::DeepTauBase
DeepTauBase(const edm::ParameterSet &cfg, const OutputCollection &outputs, const DeepTauCache *cache)
Definition: DeepTauBase.cc:76
MakerMacros.h
reco::tau::pt_weighted_dphi_strip
float pt_weighted_dphi_strip(const reco::PFTau &tau, int dm)
Definition: PFRecoTauClusterVariables.h:44
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
spr::deltaEta
static const double deltaEta
Definition: CaloConstants.h:8
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
reco::tau::countMatches
void countMatches(const reco::Muon &muon, std::vector< int > &numMatchesDT, std::vector< int > &numMatchesCSC, std::vector< int > &numMatchesRPC)
Definition: RecoTauMuonTools.cc:46
vars
vars
Definition: DeepTauId.cc:158
PVValHelper::eta
Definition: PVValidationHelpers.h:69
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
DeepTauId::createConvFeatures
void createConvFeatures(const TauType &tau, const reco::Vertex &pv, double rho, const pat::ElectronCollection &electrons, const pat::MuonCollection &muons, const pat::PackedCandidateCollection &pfCands, const CellGrid &grid, bool is_inner)
Definition: DeepTauId.cc:1100
HLT_2018_cff.dPhi
dPhi
Definition: HLT_2018_cff.py:12290
dnn_inputs_2017_v2::HadronBlockInputs
Definition: DeepTauId.cc:372
pat::MuonCollection
std::vector< Muon > MuonCollection
Definition: Muon.h:35
trackingPlots.other
other
Definition: trackingPlots.py:1465
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
DeepTauId::getValueLinear
static float getValueLinear(T value, float min_value, float max_value, bool positive)
Definition: DeepTauId.cc:907
MuonSubdetId::DT
static constexpr int DT
Definition: MuonSubdetId.h:11
DeepTauId::calculateGottfriedJacksonAngleDifference
static bool calculateGottfriedJacksonAngleDifference(const pat::Tau &tau, double &gj_diff)
Definition: DeepTauId.cc:2039
dqmdumpme.k
k
Definition: dqmdumpme.py:60
hcaldqm::quantity::min_value
const std::map< ValueQuantityType, double > min_value
Definition: ValueQuantity.h:130
CommPDSkim_cfg.maxDeltaEta
maxDeltaEta
Definition: CommPDSkim_cfg.py:75
DeepTauId::convTensor_
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > convTensor_
Definition: DeepTauId.cc:2091
DeepTauId::calculateElectronClusterVarsV2
static bool calculateElectronClusterVarsV2(const pat::Electron &ele, float &cc_ele_energy, float &cc_gamma_energy, int &cc_n_gamma)
Definition: DeepTauId.cc:923
DeepTauId::getPredictionsV2
void getPredictionsV2(const TauType &tau, const pat::ElectronCollection &electrons, const pat::MuonCollection &muons, const pat::PackedCandidateCollection &pfCands, const reco::Vertex &pv, double rho, std::vector< tensorflow::Tensor > &pred_vector)
Definition: DeepTauId.cc:1014
taus_cff.decayMode
decayMode
Definition: taus_cff.py:60
utilities.cache
def cache(function)
Definition: utilities.py:3
getGTfromDQMFile.obj
obj
Definition: getGTfromDQMFile.py:32
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
runTauDisplay.tau_phi
tau_phi
Definition: runTauDisplay.py:127
DeepTauId::version
const unsigned version
Definition: DeepTauId.cc:2087
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
DDAxes::rho
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
deep_tau::DeepTauBase::MuonCollection
pat::MuonCollection MuonCollection
Definition: DeepTauBase.h:66
TrackingRecHit::bad
Definition: TrackingRecHit.h:49
DeepTauId::fillGrids
void fillGrids(const TauType &tau, const Collection &objects, CellGrid &inner_grid, CellGrid &outer_grid)
Definition: DeepTauId.cc:1040
deep_tau::DeepTauBase::outputs_
OutputCollection outputs_
Definition: DeepTauBase.h:103
DeepTauId::getPredictions
tensorflow::Tensor getPredictions(edm::Event &event, const edm::EventSetup &es, edm::Handle< TauCollection > taus) override
Definition: DeepTauId.cc:968
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
DeepTauId::output_layer_
std::string output_layer_
Definition: DeepTauId.cc:2086
edm::ParameterSet
Definition: ParameterSet.h:36
Cell
Definition: LogEleMapdb.h:8
DeepTauId::calculateGottfriedJacksonAngleDifference
static float calculateGottfriedJacksonAngleDifference(const pat::Tau &tau)
Definition: DeepTauId.cc:2057
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
KineDebug3::count
void count()
Definition: KinematicConstrainedVertexUpdatorT.h:21
pat::PackedCandidate
Definition: PackedCandidate.h:22
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
dnn_inputs_2017_v2::EgammaBlockInputs
Definition: DeepTauId.cc:210
DeepTauId
Definition: DeepTauId.cc:808
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
reco::tau::n_photons_total
unsigned int n_photons_total(const reco::PFTau &tau)
return total number of pf photon candidates with pT>500 MeV, which are associated to signal
Definition: PFRecoTauClusterVariables.cc:157
cand
Definition: decayParser.h:34
electrons_cff.ip3d
ip3d
Definition: electrons_cff.py:363
runTauDisplay.tau_mass
tau_mass
Definition: runTauDisplay.py:128
createfilelist.int
int
Definition: createfilelist.py:10
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
dumpRecoGeometry_cfg.Muon
Muon
Definition: dumpRecoGeometry_cfg.py:190
DeepTauId::getPredictionsV1
void getPredictionsV1(const TauType &tau, const pat::ElectronCollection &electrons, const pat::MuonCollection &muons, std::vector< tensorflow::Tensor > &pred_vector)
Definition: DeepTauId.cc:1006
p4
double p4[4]
Definition: TauolaWrapper.h:92
value
Definition: value.py:1
pat::Muon::PV2D
Definition: Muon.h:235
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:50
PixelMapPlotter.inputs
inputs
Definition: PixelMapPlotter.py:490
DeepTauId::GetOutputs
static const OutputCollection & GetOutputs()
Definition: DeepTauId.cc:812
DeepTauId::createMuonBlockInputs
void createMuonBlockInputs(unsigned idx, const TauType &tau, const reco::Vertex &pv, double rho, const pat::MuonCollection &muons, const pat::PackedCandidateCollection &pfCands, const Cell &cell_map, bool is_inner)
Definition: DeepTauId.cc:1489
edm::EventSetup
Definition: EventSetup.h:57
DeepTauId::checkInputs
void checkInputs(const tensorflow::Tensor &inputs, const char *block_name, int n_inputs, int n_eta=1, int n_phi=1) const
Definition: DeepTauId.cc:946
DeepTauId::globalEndJob
static void globalEndJob(const deep_tau::DeepTauCache *cache_)
Definition: DeepTauId.cc:896
beam_dqm_sourceclient-live_cfg.minPt
minPt
Definition: beam_dqm_sourceclient-live_cfg.py:315
DeepTauId::eGammaTensor_
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > eGammaTensor_
Definition: DeepTauId.cc:2091
get
#define get
DeepTauId::tauBlockTensor_
std::unique_ptr< tensorflow::Tensor > tauBlockTensor_
Definition: DeepTauId.cc:2090
DeepTauId::muons_token_
edm::EDGetTokenT< MuonCollection > muons_token_
Definition: DeepTauId.cc:2084
TrackingRecHit::valid
Definition: TrackingRecHit.h:46
reco::HitPattern::TRACK_HITS
Definition: HitPattern.h:155
looper.cfg
cfg
Definition: looper.py:297
pat::PackedCandidateCollection
std::vector< pat::PackedCandidate > PackedCandidateCollection
Definition: PackedCandidate.h:1130
deep_tau::DeepTauBase::ElectronCollection
pat::ElectronCollection ElectronCollection
Definition: DeepTauBase.h:65
nanoDQM_cff.Electron
Electron
Definition: nanoDQM_cff.py:62
DDAxes::phi
hcaldqm::quantity::max_value
const std::map< ValueQuantityType, double > max_value
Definition: ValueQuantity.h:190
DeepTauId::electrons_token_
edm::EDGetTokenT< ElectronCollection > electrons_token_
Definition: DeepTauId.cc:2083
PVValHelper::dxy
Definition: PVValidationHelpers.h:47
pat::Electron::superCluster
reco::SuperClusterRef superCluster() const override
override the reco::GsfElectron::superCluster method, to access the internal storage of the superclust...
isFinite.h
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
DeepTauId::muonTensor_
std::array< std::unique_ptr< tensorflow::Tensor >, 2 > muonTensor_
Definition: DeepTauId.cc:2091
PVValHelper::dz
Definition: PVValidationHelpers.h:50
MuonSubdetId::RPC
static constexpr int RPC
Definition: MuonSubdetId.h:13
T
long double T
Definition: Basic3DVectorLD.h:48
HLT_2018_cff.candidates
candidates
Definition: HLT_2018_cff.py:53513
CellObjectType
CellObjectType
Definition: DeepTauId.cc:703
reco::CandidateCollection
edm::OwnVector< Candidate > CandidateCollection
collection of Candidate objects
Definition: CandidateFwd.h:21
runTauDisplay.tau_pt
tau_pt
Definition: runTauDisplay.py:125
Exception
Definition: hltDiff.cc:246
postprocess-scan-build.cells
cells
Definition: postprocess-scan-build.py:13
reco::tau::pt_weighted_dr_signal
float pt_weighted_dr_signal(const reco::PFTau &tau, int dm)
Definition: PFRecoTauClusterVariables.h:32
tensorflow::run
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:211
EgHLTOffHistBins_cfi.mass
mass
Definition: EgHLTOffHistBins_cfi.py:34
DeepTauId::createHadronsBlockInputs
void createHadronsBlockInputs(unsigned idx, const TauType &tau, const reco::Vertex &pv, double rho, const pat::PackedCandidateCollection &pfCands, const Cell &cell_map, bool is_inner)
Definition: DeepTauId.cc:1629
CommPDSkim_cfg.maxDeltaPhi
maxDeltaPhi
Definition: CommPDSkim_cfg.py:74
DeepTauId::input_layer_
std::string input_layer_
Definition: DeepTauId.cc:2086
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
reco::tau::countHits
void countHits(const reco::Muon &muon, std::vector< int > &numHitsDT, std::vector< int > &numHitsCSC, std::vector< int > &numHitsRPC)
Definition: RecoTauMuonTools.cc:7
reco::tau::pt_weighted_dr_iso
float pt_weighted_dr_iso(const reco::PFTau &tau, int dm)
Definition: PFRecoTauClusterVariables.h:52
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
DeepTauId::createEgammaBlockInputs
void createEgammaBlockInputs(unsigned idx, const TauType &tau, const reco::Vertex &pv, double rho, const pat::ElectronCollection &electrons, const pat::PackedCandidateCollection &pfCands, const Cell &cell_map, bool is_inner)
Definition: DeepTauId.cc:1270
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
cms::Exception
Definition: Exception.h:70
pat::Electron
Analysis-level electron class.
Definition: Electron.h:51
deep_tau::DeepTauCache
Definition: DeepTauBase.h:40
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
operator<
bool operator<(DTCELinkId const &lhs, DTCELinkId const &rhs)
Definition: DTCELinkId.h:70
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition: HGC3DClusterGenMatchSelector_cfi.py:7
operator[]
T operator[](int i) const
Definition: extBasic3DVector.h:222
DeepTauId::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: DeepTauId.cc:822
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
dnn_inputs_2017_v2::TauBlockInputs
Definition: DeepTauId.cc:157
DeepTauId::processIsolationPFComponents
static void processIsolationPFComponents(const pat::Tau &tau, const CandidateCollection &candidates, LorentzVectorXYZ &p4, float &pt, float &d_eta, float &d_phi, float &m, float &n)
Definition: DeepTauId.cc:2010
edm::InputTag
Definition: InputTag.h:15
begin
#define begin
Definition: vmac.h:32
DeepTauId::getPartialPredictions
tensorflow::Tensor getPartialPredictions(bool is_inner)
Definition: DeepTauId.cc:1076
deep_tau::DeepTauCache::getGraph
const tensorflow::GraphDef & getGraph(const std::string &name="") const
Definition: DeepTauBase.h:50
Output
#define Output(cl)
Definition: vmac.h:194
reco::Vertex
Definition: Vertex.h:35
hit
Definition: SiStripHitEffFromCalibTree.cc:88
HGVHistoProducerAlgoBlock_cfi.maxX
maxX
Definition: HGVHistoProducerAlgoBlock_cfi.py:154
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
pwdgSkimBPark_cfi.vertices
vertices
Definition: pwdgSkimBPark_cfi.py:7