CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
lowptgsfeleid::HeavyObjectCache Class Reference

#include <LowPtGsfElectronIDHeavyObjectCache.h>

Public Member Functions

double eval (const std::string &name, const reco::GsfElectronRef &, double rho) const
 
 HeavyObjectCache (const edm::ParameterSet &)
 
std::vector< std::string > modelNames () const
 

Private Attributes

std::vector< std::unique_ptr< const GBRForest > > models_
 
std::vector< std::string > names_
 
std::vector< double > thresholds_
 

Detailed Description

Definition at line 52 of file LowPtGsfElectronIDHeavyObjectCache.h.

Constructor & Destructor Documentation

lowptgsfeleid::HeavyObjectCache::HeavyObjectCache ( const edm::ParameterSet conf)

Definition at line 102 of file LowPtGsfElectronIDHeavyObjectCache.cc.

References createGBRForest(), Exception, edm::ParameterSet::getParameter(), dataset::name, scrapingFilter_cfi::thresh, and HGCalRecHit_cfi::weights.

102  {
103  for ( auto& name : conf.getParameter< std::vector<std::string> >("ModelNames") )
104  {
105  names_.push_back(name);
106  }
107  for ( auto& weights : conf.getParameter< std::vector<std::string> >("ModelWeights") )
108  {
110  }
111  for ( auto& thresh : conf.getParameter< std::vector<double> >("ModelThresholds") )
112  {
113  thresholds_.push_back(thresh);
114  }
115  if ( names_.size() != models_.size() ) {
116  throw cms::Exception("Incorrect configuration")
117  << "'ModelNames' size (" << names_.size()
118  << ") != 'ModelWeights' size (" << models_.size()
119  << ").\n";
120  }
121  if ( models_.size() != thresholds_.size() ) {
122  throw cms::Exception("Incorrect configuration")
123  << "'ModelWeights' size (" << models_.size()
124  << ") != 'ModelThresholds' size (" << thresholds_.size()
125  << ").\n";
126  }
127 
128  }
T getParameter(std::string const &) const
std::vector< std::unique_ptr< const GBRForest > > models_
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)

Member Function Documentation

double lowptgsfeleid::HeavyObjectCache::eval ( const std::string &  name,
const reco::GsfElectronRef ele,
double  rho 
) const

Definition at line 132 of file LowPtGsfElectronIDHeavyObjectCache.cc.

References SoftLeptonByDistance_cfi::distance, Exception, postprocess-scan-build::features, spr::find(), haddnano::inputs, and dataset::name.

135  {
136  std::vector<std::string>::const_iterator iter = std::find( names_.begin(),
137  names_.end(),
138  name );
139  if ( iter != names_.end() ) {
140  int index = std::distance(names_.begin(),iter);
141  Features features;
142  features.set(ele,rho);
143  std::vector<float> inputs = features.get();
144  return models_.at(index)->GetResponse( inputs.data() );
145  } else {
146  throw cms::Exception("Unknown model name")
147  << "'Name given: '" << name
148  << "'. Check against configuration file.\n";
149  }
150  return 0.;
151  }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< std::unique_ptr< const GBRForest > > models_
std::vector<std::string> lowptgsfeleid::HeavyObjectCache::modelNames ( ) const
inline

Member Data Documentation

std::vector< std::unique_ptr<const GBRForest> > lowptgsfeleid::HeavyObjectCache::models_
private

Definition at line 65 of file LowPtGsfElectronIDHeavyObjectCache.h.

std::vector<std::string> lowptgsfeleid::HeavyObjectCache::names_
private

Definition at line 64 of file LowPtGsfElectronIDHeavyObjectCache.h.

std::vector<double> lowptgsfeleid::HeavyObjectCache::thresholds_
private

Definition at line 66 of file LowPtGsfElectronIDHeavyObjectCache.h.