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 
9 // system include files
10 #include <memory>
11 
12 // user include files
15 
18 
25 //
26 // class declaration
27 //
28 
30 public:
31  explicit ElectronIdMVABased(const edm::ParameterSet&);
32  ~ElectronIdMVABased() override {}
33 
34 private:
35  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
36 
37  // ----------member data ---------------------------
40  const std::vector<std::string> mvaWeightFileEleID;
42  const double thresholdBarrel;
43  const double thresholdEndcap;
44  const double thresholdIsoBarrel;
45  const double thresholdIsoEndcap;
46 
47  const std::unique_ptr<const ElectronMVAEstimator> mvaID_;
48 };
49 
50 // constructor
52  : vertexToken(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexTag")))
53  , electronToken(consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("electronTag")))
54  , thresholdBarrel (iConfig.getParameter<double>("thresholdBarrel"))
55  , thresholdEndcap (iConfig.getParameter<double>("thresholdEndcap"))
56  , thresholdIsoBarrel(iConfig.getParameter<double>("thresholdIsoDR03Barrel"))
57  , thresholdIsoEndcap(iConfig.getParameter<double>("thresholdIsoDR03Endcap"))
59  .vweightsfiles = iConfig.getParameter<std::vector<std::string> >("HZZmvaWeightFile")}))
60 {
61  produces<reco::GsfElectronCollection>();
62 }
63 
64 // ------------ method called on each new Event ------------
66 {
67 
68  constexpr double etaEBEE = 1.485;
69 
70  auto mvaElectrons = std::make_unique<reco::GsfElectronCollection>();
71 
73  iEvent.getByToken(vertexToken, vertexCollection);
74  int nVtx = vertexCollection->size();
75 
77  iEvent.getByToken(electronToken,egCollection);
78  const reco::GsfElectronCollection egCandidates = (*egCollection.product());
79  for ( reco::GsfElectronCollection::const_iterator egIter = egCandidates.begin(); egIter != egCandidates.end(); ++egIter) {
80  double mvaVal = mvaID_->mva( *egIter, nVtx );
81  double isoDr03 = egIter->dr03TkSumPt() + egIter->dr03EcalRecHitSumEt() + egIter->dr03HcalTowerSumEt();
82  double eleEta = fabs(egIter->eta());
83  if (eleEta <= etaEBEE && mvaVal > thresholdBarrel && isoDr03 < thresholdIsoBarrel) {
84  mvaElectrons->push_back( *egIter );
85  reco::GsfElectron::MvaOutput myMvaOutput;
86  myMvaOutput.mva_Isolated = mvaVal;
87  mvaElectrons->back().setMvaOutput(myMvaOutput);
88  }
89  else if (eleEta > etaEBEE && mvaVal > thresholdEndcap && isoDr03 < thresholdIsoEndcap) {
90  mvaElectrons->push_back( *egIter );
91  reco::GsfElectron::MvaOutput myMvaOutput;
92  myMvaOutput.mva_Isolated = mvaVal;
93  mvaElectrons->back().setMvaOutput(myMvaOutput);
94  }
95  }
96 
97  iEvent.put(std::move(mvaElectrons));
98 
99 }
100 
101 //define this as a plug-in
const std::string path_mvaWeightFileEleID
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
const std::vector< std::string > mvaWeightFileEleID
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const std::unique_ptr< const ElectronMVAEstimator > mvaID_
const double thresholdIsoBarrel
const double thresholdIsoEndcap
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EDGetTokenT< reco::VertexCollection > vertexToken
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const double thresholdBarrel
const double thresholdEndcap
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
T const * product() const
Definition: Handle.h:74
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
#define constexpr