CMS 3D CMS Logo

LowPtGsfElectronIDProducer.cc
Go to the documentation of this file.
24 #include <string>
25 #include <vector>
26 
28 //
30 public:
32 
33  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
34 
36 
37 private:
38  double eval(const std::string& name, const reco::GsfElectronRef&, double rho, float unbiased, float field_z) const;
39 
43  const std::vector<std::string> names_;
44  const bool passThrough_;
45  const double minPtThreshold_;
46  const double maxPtThreshold_;
47  std::vector<std::unique_ptr<const GBRForest> > models_;
48  const std::vector<double> thresholds_;
50 };
51 
53 //
55  : electrons_(consumes<reco::GsfElectronCollection>(conf.getParameter<edm::InputTag>("electrons"))),
56  rho_(consumes<double>(conf.getParameter<edm::InputTag>("rho"))),
57  unbiased_(consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("unbiased"))),
58  names_(conf.getParameter<std::vector<std::string> >("ModelNames")),
59  passThrough_(conf.getParameter<bool>("PassThrough")),
60  minPtThreshold_(conf.getParameter<double>("MinPtThreshold")),
61  maxPtThreshold_(conf.getParameter<double>("MaxPtThreshold")),
62  thresholds_(conf.getParameter<std::vector<double> >("ModelThresholds")),
63  version_(conf.getParameter<std::string>("Version")) {
64  for (auto& weights : conf.getParameter<std::vector<std::string> >("ModelWeights")) {
66  }
67  if (names_.size() != models_.size()) {
68  throw cms::Exception("Incorrect configuration")
69  << "'ModelNames' size (" << names_.size() << ") != 'ModelWeights' size (" << models_.size() << ").\n";
70  }
71  if (models_.size() != thresholds_.size()) {
72  throw cms::Exception("Incorrect configuration")
73  << "'ModelWeights' size (" << models_.size() << ") != 'ModelThresholds' size (" << thresholds_.size() << ").\n";
74  }
75  if (version_ != "V0" && version_ != "V1" && !version_.empty()) {
76  throw cms::Exception("Incorrect configuration") << "Unknown Version: " << version_ << "\n";
77  }
78  for (const auto& name : names_) {
79  produces<edm::ValueMap<float> >(name);
80  }
81 }
82 
84 //
86  // Get z-component of B field
88  setup.get<IdealMagneticFieldRecord>().get(field);
89  math::XYZVector zfield(field->inTesla(GlobalPoint(0, 0, 0)));
90 
91  // Pileup
93  event.getByToken(rho_, rho);
94  if (!rho.isValid()) {
95  std::ostringstream os;
96  os << "Problem accessing rho collection for low-pT electrons" << std::endl;
97  throw cms::Exception("InvalidHandle", os.str());
98  }
99 
100  // Retrieve GsfElectrons from Event
102  event.getByToken(electrons_, electrons);
103  if (!electrons.isValid()) {
104  std::ostringstream os;
105  os << "Problem accessing low-pT electrons collection" << std::endl;
106  throw cms::Exception("InvalidHandle", os.str());
107  }
108 
109  // ElectronSeed unbiased BDT
111  event.getByToken(unbiased_, unbiasedH);
112 
113  // Iterate through Electrons, evaluate BDT, and store result
114  std::vector<std::vector<float> > output;
115  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
116  output.emplace_back(electrons->size(), -999.);
117  }
118  for (unsigned int iele = 0; iele < electrons->size(); iele++) {
119  reco::GsfElectronRef ele(electrons, iele);
120 
121  if (ele->core().isNull()) {
122  continue;
123  }
124  const auto& gsf = ele->core()->gsfTrack(); // reco::GsfTrackRef
125  if (gsf.isNull()) {
126  continue;
127  }
128  float unbiased = (*unbiasedH)[gsf];
129 
130  //if ( !passThrough_ && ( ele->pt() < minPtThreshold_ ) ) { continue; }
131  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
132  output[iname][iele] = eval(names_[iname], ele, *rho, unbiased, zfield.z());
133  }
134  }
135 
136  // Create and put ValueMap in Event
137  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
138  auto ptr = std::make_unique<edm::ValueMap<float> >(edm::ValueMap<float>());
140  filler.insert(electrons, output[iname].begin(), output[iname].end());
141  filler.fill();
142  event.put(std::move(ptr), names_[iname]);
143  }
144 }
145 
147 //
149  const std::string& name, const reco::GsfElectronRef& ele, double rho, float unbiased, float field_z) const {
150  auto iter = std::find(names_.begin(), names_.end(), name);
151  if (iter != names_.end()) {
152  int index = std::distance(names_.begin(), iter);
153  std::vector<float> inputs;
154  if (version_.empty()) { // Original XML model
156  features.set(ele, rho);
157  inputs = features.get();
158  } else if (version_ == "V0") {
159  inputs = lowptgsfeleid::features_V0(*ele, rho, unbiased);
160  } else if (version_ == "V1") {
161  inputs = lowptgsfeleid::features_V1(*ele, rho, unbiased, field_z);
162  }
163  return models_.at(index)->GetResponse(inputs.data());
164  } else {
165  throw cms::Exception("Unknown model name") << "'Name given: '" << name << "'. Check against configuration file.\n";
166  }
167  return 0.;
168 }
169 
171 //
174  desc.add<edm::InputTag>("electrons", edm::InputTag("lowPtGsfElectrons"));
175  desc.add<edm::InputTag>("unbiased", edm::InputTag("lowPtGsfElectronSeedValueMaps:unbiased"));
176  desc.add<edm::InputTag>("rho", edm::InputTag("fixedGridRhoFastjetAllTmp"));
177  desc.add<std::vector<std::string> >("ModelNames", std::vector<std::string>());
178  desc.add<std::vector<std::string> >("ModelWeights", std::vector<std::string>());
179  desc.add<std::vector<double> >("ModelThresholds", std::vector<double>());
180  desc.add<bool>("PassThrough", false);
181  desc.add<double>("MinPtThreshold", 0.5);
182  desc.add<double>("MaxPtThreshold", 15.);
183  desc.add<std::string>("Version", "");
184  descriptions.add("defaultLowPtGsfElectronID", desc);
185 }
186 
188 //
const edm::EDGetTokenT< edm::ValueMap< float > > unbiased_
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &)
double eval(const std::string &name, const reco::GsfElectronRef &, double rho, float unbiased, float field_z) const
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
const std::vector< double > thresholds_
const std::vector< std::string > names_
std::vector< std::unique_ptr< const GBRForest > > models_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
const edm::EDGetTokenT< reco::GsfElectronCollection > electrons_
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
LowPtGsfElectronIDProducer(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
#define end
Definition: vmac.h:39
std::vector< float > features_V0(reco::GsfElectron const &ele, float rho, float unbiased)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isNull() const
Checks for null.
Definition: Ref.h:248
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
const edm::EDGetTokenT< double > rho_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
T get() const
Definition: EventSetup.h:71
void set(const reco::GsfElectronRef &ele, double rho)
std::vector< float > features_V1(reco::GsfElectron const &ele, float rho, float unbiased, float field_z)
def move(src, dest)
Definition: eostools.py:511
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)
Definition: event.py:1