CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Classes | Functions
lowptgsfeleseed Namespace Reference

Classes

class  HeavyObjectCache
 

Functions

std::vector< float > features (const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
 

Function Documentation

std::vector< float > lowptgsfeleseed::features ( const reco::PreId ecal,
const reco::PreId hcal,
double  rho,
const reco::BeamSpot spot,
noZS::EcalClusterLazyTools ecalTools 
)

Definition at line 17 of file LowPtGsfElectronFeatures.cc.

References reco::PreId::chi2Ratio(), reco::PreId::clusterRef(), reco::PreId::dpt(), EcalClusterLazyToolsT< ClusterTools >::e1x5(), EcalClusterLazyToolsT< ClusterTools >::e3x3(), EcalClusterLazyToolsT< ClusterTools >::e5x5(), reco::PreId::geomMatching(), reco::PreId::gsfChi2(), EcalClusterLazyToolsT< ClusterTools >::localCovariances(), convertSQLitetoXML_cfg::output, reco::TrackBase::qualityByName(), and reco::PreId::trackRef().

Referenced by ticl::PatternRecognitionbyCLUE3D< TILES >::energyRegressionAndID(), ticl::PatternRecognitionbyCA< TILES >::energyRegressionAndID(), TrackstersMergeProducer::energyRegressionAndID(), lowptgsfeleseed::HeavyObjectCache::eval(), btagbtvdeep::getDRSubjetFeatures(), egammaTools::EgammaDNNHelper::initScalerFiles(), DeepFlavourONNXJetTagsProducer::make_inputs(), DeepVertexONNXJetTagsProducer::make_inputs(), DeepCombinedONNXJetTagsProducer::make_inputs(), DeepDoubleXONNXJetTagsProducer::make_inputs(), DeepBoostedJetTagInfoProducer::produce(), DeepDoubleXTagInfoProducer::produce(), and DeepFlavourTagInfoProducer::produce().

21  {
22  float trk_pt_ = -1.;
23  float trk_eta_ = -1.;
24  float trk_phi_ = -1.;
25  float trk_p_ = -1.;
26  float trk_nhits_ = -1.;
27  float trk_high_quality_ = -1.;
28  float trk_chi2red_ = -1.;
29  float rho_ = -1.;
30  float ktf_ecal_cluster_e_ = -1.;
31  float ktf_ecal_cluster_deta_ = -42.;
32  float ktf_ecal_cluster_dphi_ = -42.;
33  float ktf_ecal_cluster_e3x3_ = -1.;
34  float ktf_ecal_cluster_e5x5_ = -1.;
35  float ktf_ecal_cluster_covEtaEta_ = -42.;
36  float ktf_ecal_cluster_covEtaPhi_ = -42.;
37  float ktf_ecal_cluster_covPhiPhi_ = -42.;
38  float ktf_ecal_cluster_r9_ = -0.1;
39  float ktf_ecal_cluster_circularity_ = -0.1;
40  float ktf_hcal_cluster_e_ = -1.;
41  float ktf_hcal_cluster_deta_ = -42.;
42  float ktf_hcal_cluster_dphi_ = -42.;
43  float preid_gsf_dpt_ = -1.;
44  float preid_trk_gsf_chiratio_ = -1.;
45  float preid_gsf_chi2red_ = -1.;
46  float trk_dxy_sig_ = -1.; // must be last (not used by unbiased model)
47 
48  // Tracks
49  const auto& trk = ecal.trackRef(); // reco::TrackRef
50  if (trk.isNonnull()) {
51  trk_pt_ = trk->pt();
52  trk_eta_ = trk->eta();
53  trk_phi_ = trk->phi();
54  trk_p_ = trk->p();
55  trk_nhits_ = static_cast<float>(trk->found());
56  trk_high_quality_ = static_cast<float>(trk->quality(reco::TrackBase::qualityByName("highPurity")));
57  trk_chi2red_ = trk->normalizedChi2();
58  if (trk->dxy(spot) > 0.) {
59  trk_dxy_sig_ = trk->dxyError() / trk->dxy(spot); //@@ to be consistent with the training based on 94X MC
60  }
61  ktf_ecal_cluster_dphi_ *= trk->charge(); //@@ to be consistent with the training based on 94X MC
62  }
63 
64  // Rho
65  rho_ = static_cast<float>(rho);
66 
67  // ECAL clusters
68  const auto& ecal_clu = ecal.clusterRef(); // reco::PFClusterRef
69  if (ecal_clu.isNonnull()) {
70  ktf_ecal_cluster_e_ = ecal_clu->energy();
71  ktf_ecal_cluster_deta_ = ecal.geomMatching()[0];
72  ktf_ecal_cluster_dphi_ = ecal.geomMatching()[1];
73  ktf_ecal_cluster_e3x3_ = tools.e3x3(*ecal_clu);
74  ktf_ecal_cluster_e5x5_ = tools.e5x5(*ecal_clu);
75  const auto& covs = tools.localCovariances(*ecal_clu);
76  ktf_ecal_cluster_covEtaEta_ = covs[0];
77  ktf_ecal_cluster_covEtaPhi_ = covs[1];
78  ktf_ecal_cluster_covPhiPhi_ = covs[2];
79  if (ktf_ecal_cluster_e_ > 0.) {
80  ktf_ecal_cluster_r9_ = ktf_ecal_cluster_e3x3_ / ktf_ecal_cluster_e_;
81  }
82  if (ktf_ecal_cluster_e5x5_ > 0.) {
83  ktf_ecal_cluster_circularity_ = 1. - tools.e1x5(*ecal_clu) / ktf_ecal_cluster_e5x5_;
84  } else {
85  ktf_ecal_cluster_circularity_ = -0.1;
86  }
87  }
88 
89  // HCAL clusters
90  const auto& hcal_clu = hcal.clusterRef(); // reco::PFClusterRef
91  if (hcal_clu.isNonnull()) {
92  ktf_hcal_cluster_e_ = hcal_clu->energy();
93  ktf_hcal_cluster_deta_ = hcal.geomMatching()[0];
94  ktf_hcal_cluster_dphi_ = hcal.geomMatching()[1];
95  }
96 
97  // PreId
98  preid_gsf_dpt_ = ecal.dpt();
99  preid_trk_gsf_chiratio_ = ecal.chi2Ratio();
100  preid_gsf_chi2red_ = ecal.gsfChi2();
101 
102  // Set contents of vector
103  std::vector<float> output = {trk_pt_,
104  trk_eta_,
105  trk_phi_,
106  trk_p_,
107  trk_nhits_,
108  trk_high_quality_,
109  trk_chi2red_,
110  rho_,
111  ktf_ecal_cluster_e_,
112  ktf_ecal_cluster_deta_,
113  ktf_ecal_cluster_dphi_,
114  ktf_ecal_cluster_e3x3_,
115  ktf_ecal_cluster_e5x5_,
116  ktf_ecal_cluster_covEtaEta_,
117  ktf_ecal_cluster_covEtaPhi_,
118  ktf_ecal_cluster_covPhiPhi_,
119  ktf_ecal_cluster_r9_,
120  ktf_ecal_cluster_circularity_,
121  ktf_hcal_cluster_e_,
122  ktf_hcal_cluster_deta_,
123  ktf_hcal_cluster_dphi_,
124  preid_gsf_dpt_,
125  preid_trk_gsf_chiratio_,
126  preid_gsf_chi2red_,
127  trk_dxy_sig_};
128  return output;
129  };
reco::TrackRef trackRef() const
Definition: PreId.h:99
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
PFClusterRef clusterRef() const
Definition: PreId.h:100
float gsfChi2() const
Definition: PreId.h:89
const std::vector< float > & geomMatching() const
Access methods.
Definition: PreId.h:78
float dpt() const
Definition: PreId.h:98
float chi2Ratio() const
Definition: PreId.h:88