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 edm::Ptr<reco::GsfElectron>&, double rho, float unbiased, float field_z) const;
39 
40  const bool usePAT_;
45  const std::vector<std::string> names_;
46  const bool passThrough_;
47  const double minPtThreshold_;
48  const double maxPtThreshold_;
49  std::vector<std::unique_ptr<const GBRForest> > models_;
50  const std::vector<double> thresholds_;
52 };
53 
55 //
57  : usePAT_(conf.getParameter<bool>("usePAT")),
58  electrons_(),
59  patElectrons_(),
60  rho_(consumes<double>(conf.getParameter<edm::InputTag>("rho"))),
61  unbiased_(),
62  names_(conf.getParameter<std::vector<std::string> >("ModelNames")),
63  passThrough_(conf.getParameter<bool>("PassThrough")),
64  minPtThreshold_(conf.getParameter<double>("MinPtThreshold")),
65  maxPtThreshold_(conf.getParameter<double>("MaxPtThreshold")),
66  thresholds_(conf.getParameter<std::vector<double> >("ModelThresholds")),
67  version_(conf.getParameter<std::string>("Version")) {
68  if (usePAT_) {
69  patElectrons_ = consumes<pat::ElectronCollection>(conf.getParameter<edm::InputTag>("electrons"));
70  } else {
71  electrons_ = consumes<reco::GsfElectronCollection>(conf.getParameter<edm::InputTag>("electrons"));
72  unbiased_ = consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("unbiased"));
73  }
74  for (auto& weights : conf.getParameter<std::vector<std::string> >("ModelWeights")) {
76  }
77  if (names_.size() != models_.size()) {
78  throw cms::Exception("Incorrect configuration")
79  << "'ModelNames' size (" << names_.size() << ") != 'ModelWeights' size (" << models_.size() << ").\n";
80  }
81  if (models_.size() != thresholds_.size()) {
82  throw cms::Exception("Incorrect configuration")
83  << "'ModelWeights' size (" << models_.size() << ") != 'ModelThresholds' size (" << thresholds_.size() << ").\n";
84  }
85  if (version_ != "V0" && version_ != "V1" && !version_.empty()) {
86  throw cms::Exception("Incorrect configuration") << "Unknown Version: " << version_ << "\n";
87  }
88  for (const auto& name : names_) {
89  produces<edm::ValueMap<float> >(name);
90  }
91 }
92 
94 //
96  // Get z-component of B field
98  setup.get<IdealMagneticFieldRecord>().get(field);
99  math::XYZVector zfield(field->inTesla(GlobalPoint(0, 0, 0)));
100 
101  // Pileup
103  event.getByToken(rho_, rho);
104  if (!rho.isValid()) {
105  std::ostringstream os;
106  os << "Problem accessing rho collection for low-pT electrons" << std::endl;
107  throw cms::Exception("InvalidHandle", os.str());
108  }
109 
110  // Retrieve pat::Electrons or reco::GsfElectrons from Event
113  if (usePAT_) {
114  event.getByToken(patElectrons_, patElectrons);
115  } else {
116  event.getByToken(electrons_, electrons);
117  }
118 
119  // ElectronSeed unbiased BDT
121  if (!unbiased_.isUninitialized()) {
122  event.getByToken(unbiased_, unbiasedH);
123  }
124 
125  // Iterate through Electrons, evaluate BDT, and store result
126  std::vector<std::vector<float> > output;
127  unsigned int nElectrons = usePAT_ ? patElectrons->size() : electrons->size();
128  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
129  output.emplace_back(nElectrons, -999.);
130  }
131 
132  if (usePAT_) {
133  for (unsigned int iele = 0; iele < nElectrons; iele++) {
134  edm::Ptr<pat::Electron> ele(patElectrons, iele);
135  if (!ele->isElectronIDAvailable("unbiased")) {
136  continue;
137  }
138  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
139  output[iname][iele] = eval(names_[iname], ele, *rho, ele->electronID("unbiased"), zfield.z());
140  }
141  }
142  } else {
143  for (unsigned int iele = 0; iele < nElectrons; iele++) {
144  edm::Ptr<reco::GsfElectron> ele(electrons, iele);
145  if (ele->core().isNull()) {
146  continue;
147  }
148  const auto& gsf = ele->core()->gsfTrack(); // reco::GsfTrackRef
149  if (gsf.isNull()) {
150  continue;
151  }
152  float unbiased = (*unbiasedH)[gsf];
153  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
154  output[iname][iele] = eval(names_[iname], ele, *rho, unbiased, zfield.z());
155  }
156  }
157  }
158 
159  // Create and put ValueMap in Event
160  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
161  auto ptr = std::make_unique<edm::ValueMap<float> >(edm::ValueMap<float>());
163  if (usePAT_) {
164  filler.insert(patElectrons, output[iname].begin(), output[iname].end());
165  } else {
166  filler.insert(electrons, output[iname].begin(), output[iname].end());
167  }
168  filler.fill();
169  event.put(std::move(ptr), names_[iname]);
170  }
171 }
172 
174 //
176  const std::string& name, const edm::Ptr<reco::GsfElectron>& ele, double rho, float unbiased, float field_z) const {
177  auto iter = std::find(names_.begin(), names_.end(), name);
178  if (iter != names_.end()) {
179  int index = std::distance(names_.begin(), iter);
180  std::vector<float> inputs;
181  if (version_.empty()) { // Original XML model
183  features.set(*ele, rho);
184  inputs = features.get();
185  } else if (version_ == "V0") {
186  inputs = lowptgsfeleid::features_V0(*ele, rho, unbiased);
187  } else if (version_ == "V1") {
188  inputs = lowptgsfeleid::features_V1(*ele, rho, unbiased, field_z);
189  }
190  return models_.at(index)->GetResponse(inputs.data());
191  } else {
192  throw cms::Exception("Unknown model name") << "'Name given: '" << name << "'. Check against configuration file.\n";
193  }
194  return 0.;
195 }
196 
198 //
201  desc.add<bool>("usePAT", false);
202  desc.add<edm::InputTag>("electrons", edm::InputTag("lowPtGsfElectrons"));
203  desc.addOptional<edm::InputTag>("unbiased", edm::InputTag("lowPtGsfElectronSeedValueMaps:unbiased"));
204  desc.add<edm::InputTag>("rho", edm::InputTag("fixedGridRhoFastjetAllTmp"));
205  desc.add<std::vector<std::string> >("ModelNames", std::vector<std::string>());
206  desc.add<std::vector<std::string> >("ModelWeights", std::vector<std::string>());
207  desc.add<std::vector<double> >("ModelThresholds", std::vector<double>());
208  desc.add<bool>("PassThrough", false);
209  desc.add<double>("MinPtThreshold", 0.5);
210  desc.add<double>("MaxPtThreshold", 15.);
211  desc.add<std::string>("Version", "");
212  descriptions.add("defaultLowPtGsfElectronID", desc);
213 }
214 
216 //
T getParameter(std::string const &) const
bool isElectronIDAvailable(const std::string &name) const
Returns true if a specific ID is available in this pat::Electron.
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &)
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
edm::EDGetTokenT< reco::GsfElectronCollection > electrons_
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
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
edm::EDGetTokenT< pat::ElectronCollection > patElectrons_
float electronID(const std::string &name) const
Returns a specific electron ID associated to the pat::Electron given its name.
#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.
double eval(const std::string &name, const edm::Ptr< reco::GsfElectron > &, double rho, float unbiased, float field_z) const
const edm::EDGetTokenT< double > rho_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
edm::EDGetTokenT< edm::ValueMap< float > > unbiased_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual GsfElectronCoreRef core() const
Definition: GsfElectron.cc:8
#define begin
Definition: vmac.h:32
HLT enums.
T get() const
Definition: EventSetup.h:71
bool isUninitialized() const
Definition: EDGetToken.h:70
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