CMS 3D CMS Logo

LowPtGsfElectronSeedValueMapsProducer.cc
Go to the documentation of this file.
12 
14 //
16  : gsfTracks_(consumes<reco::GsfTrackCollection>(conf.getParameter<edm::InputTag>("gsfTracks"))),
17  preIdsValueMap_(consumes<edm::ValueMap<reco::PreIdRef> >(conf.getParameter<edm::InputTag>("preIdsValueMap"))),
18  names_(conf.getParameter<std::vector<std::string> >("ModelNames")) {
19  for (const auto& name : names_) {
20  produces<edm::ValueMap<float> >(name);
21  }
22 }
23 
25 //
27 
29 //
31  // Retrieve GsfTracks from Event
33  event.getByToken(gsfTracks_, gsfTracks);
34  if (!gsfTracks.isValid()) {
35  edm::LogError("Problem with gsfTracks handle");
36  }
37 
38  // Retrieve PreIds from Event
40  event.getByToken(preIdsValueMap_, preIdsValueMap);
41  if (!preIdsValueMap.isValid()) {
42  edm::LogError("Problem with preIdsValueMap handle");
43  }
44 
45  // Iterate through GsfTracks, extract BDT output, and store result in ValueMap for each model
46  std::vector<std::vector<float> > output;
47  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
48  output.push_back(std::vector<float>(gsfTracks->size(), -999.));
49  }
50  for (unsigned int igsf = 0; igsf < gsfTracks->size(); igsf++) {
51  reco::GsfTrackRef gsf(gsfTracks, igsf);
52  if (gsf.isNonnull() && gsf->extra().isNonnull() && gsf->extra()->seedRef().isNonnull()) {
53  reco::ElectronSeedRef seed = gsf->extra()->seedRef().castTo<reco::ElectronSeedRef>();
54  if (seed.isNonnull() && seed->ctfTrack().isNonnull()) {
55  const reco::PreIdRef preid = (*preIdsValueMap)[seed->ctfTrack()];
56  if (preid.isNonnull()) {
57  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
58  output[iname][igsf] = preid->mva(iname);
59  }
60  }
61  }
62  }
63  }
64 
65  // Create and put ValueMap in Event
66  for (unsigned int iname = 0; iname < names_.size(); ++iname) {
67  auto ptr = std::make_unique<edm::ValueMap<float> >(edm::ValueMap<float>());
69  filler.insert(gsfTracks, output[iname].begin(), output[iname].end());
70  filler.fill();
71  event.put(std::move(ptr), names_[iname]);
72  }
73 }
74 
76 //
79  desc.add<edm::InputTag>("gsfTracks", edm::InputTag("lowPtGsfEleGsfTracks"));
80  desc.add<edm::InputTag>("preIdsValueMap", edm::InputTag("lowPtGsfElectronSeeds"));
81  desc.add<std::vector<std::string> >("ModelNames", {"unbiased", "ptbiased"});
82  descriptions.add("lowPtGsfElectronSeedValueMaps", desc);
83 }
84 
86 //
const edm::EDGetTokenT< reco::GsfTrackCollection > gsfTracks_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
const edm::EDGetTokenT< edm::ValueMap< reco::PreIdRef > > preIdsValueMap_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< GsfTrack > GsfTrackCollection
collection of GsfTracks
Definition: GsfTrackFwd.h:9
#define end
Definition: vmac.h:39
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
edm::Ref< reco::PreIdCollection > PreIdRef
Definition: PreIdFwd.h:8
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &)
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
void produce(edm::Event &, const edm::EventSetup &) override