CMS 3D CMS Logo

LowPtGsfElectronIDProducer.cc
Go to the documentation of this file.
8 
10 //
13  gsfElectrons_(consumes<reco::GsfElectronCollection>(conf.getParameter<edm::InputTag>("electrons"))),
14  rho_(consumes<double>(conf.getParameter<edm::InputTag>("rho"))),
15  names_(conf.getParameter< std::vector<std::string> >("ModelNames")),
16  passThrough_(conf.getParameter<bool>("PassThrough")),
17  minPtThreshold_(conf.getParameter<double>("MinPtThreshold")),
18  maxPtThreshold_(conf.getParameter<double>("MaxPtThreshold"))
19 {
20  for ( const auto& name : names_ ) {
21  produces< edm::ValueMap<float> >(name);
22  }
23 }
24 
26 //
28 
30 //
32 
33  // Pileup
35  event.getByToken(rho_,rho);
36  if ( !rho.isValid() ) { edm::LogError("Problem with rho handle"); }
37 
38  // Retrieve GsfElectrons from Event
40  event.getByToken(gsfElectrons_,gsfElectrons);
41  if ( !gsfElectrons.isValid() ) { edm::LogError("Problem with gsfElectrons handle"); }
42 
43  // Iterate through Electrons, evaluate BDT, and store result
44  std::vector< std::vector<float> > output;
45  for ( unsigned int iname = 0; iname < names_.size(); ++iname ) {
46  output.push_back( std::vector<float>(gsfElectrons->size(),-999.) );
47  }
48  for ( unsigned int iele = 0; iele < gsfElectrons->size(); iele++ ) {
49  reco::GsfElectronRef ele(gsfElectrons,iele);
50  //if ( !passThrough_ && ( ele->pt() < minPtThreshold_ ) ) { continue; }
51  for ( unsigned int iname = 0; iname < names_.size(); ++iname ) {
52  output[iname][iele] = globalCache()->eval( names_[iname], ele, *rho );
53  }
54  }
55 
56  // Create and put ValueMap in Event
57  for ( unsigned int iname = 0; iname < names_.size(); ++iname ) {
58  auto ptr = std::make_unique< edm::ValueMap<float> >( edm::ValueMap<float>() );
60  filler.insert(gsfElectrons, output[iname].begin(), output[iname].end());
61  filler.fill();
62  reco::GsfElectronRef ele(gsfElectrons,0);
63  event.put(std::move(ptr),names_[iname]);
64  }
65 
66 }
67 
69 //
71 {
73  desc.add<edm::InputTag>("electrons",edm::InputTag("lowPtGsfElectrons"));
74  desc.add<edm::InputTag>("rho",edm::InputTag("fixedGridRhoFastjetAllTmp"));
75  desc.add< std::vector<std::string> >("ModelNames",std::vector<std::string>());
76  desc.add< std::vector<std::string> >("ModelWeights",std::vector<std::string>());
77  desc.add< std::vector<double> >("ModelThresholds",std::vector<double>());
78  desc.add<bool>("PassThrough",false);
79  desc.add<double>("MinPtThreshold",0.5);
80  desc.add<double>("MaxPtThreshold",15.);
81  descriptions.add("defaultLowPtGsfElectronID",desc);
82 }
83 
85 //
static void fillDescriptions(edm::ConfigurationDescriptions &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
LowPtGsfElectronIDProducer(const edm::ParameterSet &, const lowptgsfeleid::HeavyObjectCache *)
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
const std::vector< std::string > names_
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
const edm::EDGetTokenT< reco::GsfElectronCollection > gsfElectrons_
void produce(edm::Event &, const edm::EventSetup &) override
#define end
Definition: vmac.h:39
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
const edm::EDGetTokenT< double > rho_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1