CMS 3D CMS Logo

ElectronIdMVABased.cc
Go to the documentation of this file.
1 //
2 //
3 // Original Author: Zablocki Jakub
4 // Created: Thu Feb 9 10:47:50 CST 2012
5 //
6 //
7 
8 // system include files
9 #include <memory>
10 
11 // user include files
14 
17 
24 //
25 // class declaration
26 //
27 
29 public:
30  explicit ElectronIdMVABased(const edm::ParameterSet&);
31  ~ElectronIdMVABased() override {}
32 
33 private:
34  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
35 
36  // ----------member data ---------------------------
39  const std::vector<std::string> mvaWeightFileEleID;
41  const double thresholdBarrel;
42  const double thresholdEndcap;
43  const double thresholdIsoBarrel;
44  const double thresholdIsoEndcap;
45 
46  const std::unique_ptr<const ElectronMVAEstimator> mvaID_;
47 };
48 
49 // constructor
51  : vertexToken(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexTag"))),
52  electronToken(consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("electronTag"))),
53  thresholdBarrel(iConfig.getParameter<double>("thresholdBarrel")),
54  thresholdEndcap(iConfig.getParameter<double>("thresholdEndcap")),
55  thresholdIsoBarrel(iConfig.getParameter<double>("thresholdIsoDR03Barrel")),
56  thresholdIsoEndcap(iConfig.getParameter<double>("thresholdIsoDR03Endcap")),
58  .vweightsfiles = iConfig.getParameter<std::vector<std::string> >("HZZmvaWeightFile")})) {
59  produces<reco::GsfElectronCollection>();
60 }
61 
62 // ------------ method called on each new Event ------------
64  constexpr double etaEBEE = 1.485;
65 
66  auto mvaElectrons = std::make_unique<reco::GsfElectronCollection>();
67 
69  iEvent.getByToken(vertexToken, vertexCollection);
70  int nVtx = vertexCollection->size();
71 
73  iEvent.getByToken(electronToken, egCollection);
74  const reco::GsfElectronCollection egCandidates = (*egCollection.product());
75  for (reco::GsfElectronCollection::const_iterator egIter = egCandidates.begin(); egIter != egCandidates.end();
76  ++egIter) {
77  double mvaVal = mvaID_->mva(*egIter, nVtx);
78  double isoDr03 = egIter->dr03TkSumPt() + egIter->dr03EcalRecHitSumEt() + egIter->dr03HcalTowerSumEt();
79  double eleEta = fabs(egIter->eta());
80  if (eleEta <= etaEBEE && mvaVal > thresholdBarrel && isoDr03 < thresholdIsoBarrel) {
81  mvaElectrons->push_back(*egIter);
82  reco::GsfElectron::MvaOutput myMvaOutput;
83  myMvaOutput.mva_Isolated = mvaVal;
84  mvaElectrons->back().setMvaOutput(myMvaOutput);
85  } else if (eleEta > etaEBEE && mvaVal > thresholdEndcap && isoDr03 < thresholdIsoEndcap) {
86  mvaElectrons->push_back(*egIter);
87  reco::GsfElectron::MvaOutput myMvaOutput;
88  myMvaOutput.mva_Isolated = mvaVal;
89  mvaElectrons->back().setMvaOutput(myMvaOutput);
90  }
91  }
92 
94 }
95 
96 //define this as a plug-in
const std::string path_mvaWeightFileEleID
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const std::vector< std::string > mvaWeightFileEleID
T const * product() const
Definition: Handle.h:70
const std::unique_ptr< const ElectronMVAEstimator > mvaID_
const double thresholdIsoBarrel
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
const double thresholdIsoEndcap
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
const edm::EDGetTokenT< reco::VertexCollection > vertexToken
int iEvent
Definition: GenABIO.cc:224
const double thresholdBarrel
const double thresholdEndcap
const edm::EDGetTokenT< reco::GsfElectronCollection > electronToken
~ElectronIdMVABased() override
fixed size matrix
ElectronIdMVABased(const edm::ParameterSet &)
HLT enums.
def move(src, dest)
Definition: eostools.py:511