CMS 3D CMS Logo

LowPtGsfElectronIDHeavyObjectCache.cc
Go to the documentation of this file.
11 #include <string>
12 
13 namespace lowptgsfeleid {
14 
16  //
17  std::vector<float> Features::get() {
18  std::vector<float> output = {
19  rho_,
20  ele_pt_,
21  sc_eta_,
29  trk_nhits_,
32  brem_frac_,
33  gsf_nhits_,
39  sc_E_,
40  trk_p_,
41  };
42  return output;
43  }
44 
46  //
47  void Features::set( const reco::GsfElectronRef& ele, double rho ) {
48 
49  // KF tracks
50  if ( ele->core().isNonnull() ) {
51  reco::TrackRef trk = ele->core()->ctfTrack(); //@@ is this what we want?!
52  if ( trk.isNonnull() ) {
53  trk_p_ = float(trk->p());
54  trk_nhits_ = float(trk->found());
55  trk_chi2red_ = float(trk->normalizedChi2());
56  }
57  }
58 
59  // GSF tracks
60  if ( ele->core().isNonnull() ) {
61  reco::GsfTrackRef gsf = ele->core()->gsfTrack();
62  if ( gsf.isNonnull() ) {
63  gsf_nhits_ = gsf->found();
64  gsf_chi2red_ = gsf->normalizedChi2();
65  }
66  }
67 
68  // Super clusters
69  if ( ele->core().isNonnull() ) {
70  reco::SuperClusterRef sc = ele->core()->superCluster();
71  if ( sc.isNonnull() ) {
72  sc_E_ = sc->energy();
73  sc_eta_ = sc->eta();
74  sc_etaWidth_ = sc->etaWidth();
75  sc_phiWidth_ = sc->phiWidth();
76  }
77  }
78 
79  // Track-cluster matching
80  if ( ele.isNonnull() ) {
81  match_seed_dEta_ = ele->deltaEtaSeedClusterTrackAtCalo();
82  match_eclu_EoverP_ = (1./ele->ecalEnergy()) - (1./ele->p());
83  match_SC_EoverP_ = ele->eSuperClusterOverP();
84  match_SC_dEta_ = ele->deltaEtaSuperClusterTrackAtVtx();
85  match_SC_dPhi_ = ele->deltaPhiSuperClusterTrackAtVtx();
86  }
87 
88  // Shower shape vars
89  if ( ele.isNonnull() ) {
90  shape_full5x5_sigmaIetaIeta_ = ele->full5x5_sigmaIetaIeta();
91  shape_full5x5_sigmaIphiIphi_ = ele->full5x5_sigmaIphiIphi();
92  shape_full5x5_HoverE_ = ele->full5x5_hcalOverEcal();
93  shape_full5x5_r9_ = ele->full5x5_r9();
94  shape_full5x5_circularity_ = 1. - ele->full5x5_e1x5() / ele->full5x5_e5x5();
95  }
96 
97  // Misc
98  rho_ = rho;
99  if ( ele.isNonnull() ) {
100  brem_frac_ = ele->fbrem();
101  ele_pt_ = ele->pt();
102  }
103 
104  };
105 
107  //
109  for ( auto& name : conf.getParameter< std::vector<std::string> >("ModelNames") )
110  {
111  names_.push_back(name);
112  }
113  for ( auto& weights : conf.getParameter< std::vector<std::string> >("ModelWeights") )
114  {
115  models_.push_back(createGBRForest(edm::FileInPath(weights)));
116  }
117  for ( auto& thresh : conf.getParameter< std::vector<double> >("ModelThresholds") )
118  {
119  thresholds_.push_back(thresh);
120  }
121  if ( names_.size() != models_.size() ) {
122  throw cms::Exception("Incorrect configuration")
123  << "'ModelNames' size (" << names_.size()
124  << ") != 'ModelWeights' size (" << models_.size()
125  << ").\n";
126  }
127  if ( models_.size() != thresholds_.size() ) {
128  throw cms::Exception("Incorrect configuration")
129  << "'ModelWeights' size (" << models_.size()
130  << ") != 'ModelThresholds' size (" << thresholds_.size()
131  << ").\n";
132  }
133 
134  }
135 
137  //
139  const reco::GsfElectronRef& ele,
140  double rho ) const
141  {
142  std::vector<std::string>::const_iterator iter = std::find( names_.begin(),
143  names_.end(),
144  name );
145  if ( iter != names_.end() ) {
146  int index = std::distance(names_.begin(),iter);
148  features.set(ele,rho);
149  std::vector<float> inputs = features.get();
150  return models_.at(index)->GetResponse( inputs.data() );
151  } else {
152  throw cms::Exception("Unknown model name")
153  << "'Name given: '" << name
154  << "'. Check against configuration file.\n";
155  }
156  return 0.;
157  }
158 
159 } // namespace lowptgsfeleid
T getParameter(std::string const &) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
double eval(const std::string &name, const reco::GsfElectronRef &, double rho) const
void set(const reco::GsfElectronRef &ele, double rho)
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)