CMS 3D CMS Logo

LowPtGsfElectronSeedHeavyObjectCache.cc
Go to the documentation of this file.
10 
11 #include <string>
12 
13 namespace lowptgsfeleseed {
14 
16  //
17  std::vector<float> Features::get() {
18  std::vector<float> output = {
19  trk_pt_,
20  trk_eta_,
21  trk_phi_,
22  trk_p_,
23  trk_nhits_,
26  rho_,
44  };
45  return output;
46  }
47 
49  //
51  const reco::PreId& hcal,
52  double rho,
53  const reco::BeamSpot& spot,
55 
56  // Tracks
57  reco::TrackRef trk = ecal.trackRef();
58  if ( trk.isNonnull() ) {
59  trk_pt_ = trk->pt();
60  trk_eta_ = trk->eta();
61  trk_phi_ = trk->phi();
62  trk_p_ = trk->p();
63  trk_nhits_ = static_cast<float>(trk->found());
64  trk_high_quality_ = static_cast<float>(trk->quality(reco::TrackBase::qualityByName("highPurity")));
65  trk_chi2red_ = trk->normalizedChi2();
66  if ( trk->dxy(spot) > 0. ) {
67  trk_dxy_sig_ = trk->dxyError() / trk->dxy(spot); //@@ to be consistent with the training based on 94X MC
68  }
69  ktf_ecal_cluster_dphi_ *= trk->charge(); //@@ to be consistent with the training based on 94X MC
70  }
71 
72  // Rho
73  rho_ = static_cast<float>(rho);
74 
75  // ECAL clusters
76  reco::PFClusterRef ecal_clu = ecal.clusterRef();
77  if ( ecal_clu.isNonnull() ) {
78  ktf_ecal_cluster_e_ = ecal_clu->energy();
81  ktf_ecal_cluster_e3x3_ = tools.e3x3(*ecal_clu);
82  ktf_ecal_cluster_e5x5_ = tools.e5x5(*ecal_clu);
83  auto covs = tools.localCovariances(*ecal_clu);
87  if ( ktf_ecal_cluster_e_ > 0. ) {
89  }
90  if ( ktf_ecal_cluster_e5x5_ > 0. ) {
92  } else {
94  }
95  }
96 
97  // HCAL clusters
98  reco::PFClusterRef hcal_clu = hcal.clusterRef();
99  if ( hcal_clu.isNonnull() ) {
100  ktf_hcal_cluster_e_ = hcal_clu->energy();
103  }
104 
105  // PreId
106  preid_gsf_dpt_ = ecal.dpt();
108  preid_gsf_chi2red_ = ecal.gsfChi2();
109 
110  };
111 
113  //
115  for ( auto& name : conf.getParameter< std::vector<std::string> >("ModelNames") )
116  {
117  names_.push_back(name);
118  }
119  for ( auto& weights : conf.getParameter< std::vector<std::string> >("ModelWeights") )
120  {
121  models_.push_back(createGBRForest(edm::FileInPath(weights)));
122  }
123  for ( auto& thresh : conf.getParameter< std::vector<double> >("ModelThresholds") )
124  {
125  thresholds_.push_back(thresh);
126  }
127  if ( names_.size() != models_.size() ) {
128  throw cms::Exception("Incorrect configuration")
129  << "'ModelNames' size (" << names_.size()
130  << ") != 'ModelWeights' size (" << models_.size()
131  << ").\n";
132  }
133  if ( models_.size() != thresholds_.size() ) {
134  throw cms::Exception("Incorrect configuration")
135  << "'ModelWeights' size (" << models_.size()
136  << ") != 'ModelThresholds' size (" << thresholds_.size()
137  << ").\n";
138  }
139 
140  }
141 
143  //
145  reco::PreId& ecal,
146  reco::PreId& hcal,
147  double rho,
148  const reco::BeamSpot& spot,
149  noZS::EcalClusterLazyTools& ecalTools ) const
150  {
151  std::vector<std::string>::const_iterator iter = std::find( names_.begin(),
152  names_.end(),
153  name );
154  if ( iter != names_.end() ) {
155  int index = std::distance(names_.begin(),iter);
156  Features features;
157  features.set(ecal,hcal,rho,spot,ecalTools);
158  std::vector<float> inputs = features.get();
159  float output = models_.at(index)->GetResponse( inputs.data() );
160  bool pass = output > thresholds_.at(index);
161  ecal.setMVA(pass,output,index);
162  return pass;
163  } else {
164  throw cms::Exception("Unknown model name")
165  << "'Name given: '" << name
166  << "'. Check against configuration file.\n";
167  }
168  }
169 
170 } // namespace lowptgsfeleseed
T getParameter(std::string const &) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
float e5x5(const reco::BasicCluster &cluster)
reco::TrackRef trackRef() const
Definition: PreId.h:111
float e3x3(const reco::BasicCluster &cluster)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void setMVA(bool accepted, float mva, unsigned n=0)
Definition: PreId.h:72
void set(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
Definition: tools.py:1
bool eval(const std::string &name, reco::PreId &ecal, reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools) const
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
std::vector< float > localCovariances(const reco::BasicCluster &cluster, float w0=4.7)
PFClusterRef clusterRef() const
Definition: PreId.h:112
float gsfChi2() const
Definition: PreId.h:101
float e1x5(const reco::BasicCluster &cluster)
const std::vector< float > & geomMatching() const
Access methods.
Definition: PreId.h:90
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)
float dpt() const
Definition: PreId.h:110
float chi2Ratio() const
Definition: PreId.h:100