CMS 3D CMS Logo

EgammaProbeSelector.cc
Go to the documentation of this file.
1 
14 #include <iostream>
15 #include <TMath.h>
16 
17 using namespace edm;
18 using namespace std;
19 using namespace reco;
20 
21 
23 
24  jetToken = consumes<CaloJetCollection>(iConfig.getParameter<InputTag>("JetCollection"));
25  minNumberOfjets = iConfig.getParameter<int>("MinNumberOfJets");
26  jetEtMin = iConfig.getParameter<double>("JetEtMin");
27  jetEtaMin = iConfig.getParameter<double>("JetEtaMin");
28  jetEtaMax = iConfig.getParameter<double>("JetEtaMax");
29 
30  nEvents = 0;
31  nSelectedEvents = 0;
32 
33  scToken = consumes<SuperClusterCollection>(iConfig.getParameter<InputTag>("SuperClusterBarrelCollection"));
34  scEEToken = consumes<SuperClusterCollection>(iConfig.getParameter<InputTag>("SuperClusterEndCapCollection"));
35  minNumberOfSuperClusters = iConfig.getParameter<int>("MinNumberOfSuperClusters");
36  scEtMin = iConfig.getParameter<double>("ScEtMin");
37  scEtaMin = iConfig.getParameter<double>("ScEtaMin");
38  scEtaMax = iConfig.getParameter<double>("ScEtaMax");
39 }
40 
41 
43  edm::LogVerbatim("EgammaProbeSelector")
44  << " Number_events_read " << nEvents
45  << " Number_events_kept " << nSelectedEvents
46  << " Efficiency " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl;
47 
48 }
49 
50 
52  int nSC = 0;
53  int nJets = 0;
54  nEvents++;
55 
56  // Get reconstructed Jets, check their Et and eta, and count their number
57  Handle<CaloJetCollection> jetHandle;
58  iEvent.getByToken(jetToken,jetHandle);
59 
60  if(jetHandle.isValid()){
61  const reco::CaloJetCollection & jets = *(jetHandle.product());
62  CaloJetCollection::const_iterator ijet;
63  for(ijet = jets.begin(); ijet!= jets.end(); ijet++){
64  if(ijet->pt() > jetEtMin &&
65  ijet->eta() > jetEtaMin &&
66  ijet->eta() < jetEtaMax )
67  nJets++;
68  }
69  }
70 
71  // Get Super Clusters, check their Et and eta, and count their number
72  // EB
74  iEvent.getByToken(scToken,scHandle);
75 
76  if(scHandle.isValid()){
77  const reco::SuperClusterCollection & SCs = *(scHandle.product());
78 
79  for(reco::SuperClusterCollection::const_iterator scIt = SCs.begin();
80  scIt!= SCs.end(); scIt++){
81 
82  double scEt = scIt->energy() * TMath::Sin(scIt->position().theta());
83  double _eta = scIt->position().eta();
84  if(scEt > scEtMin && _eta > scEtaMin && _eta < scEtaMax) nSC++;
85  }
86  }
87 
88  // EE
90  iEvent.getByToken(scEEToken,scEEHandle);
91 
92  if(scEEHandle.isValid()){
93  const reco::SuperClusterCollection & SCs = *(scEEHandle.product());
94  for(reco::SuperClusterCollection::const_iterator scIt = SCs.begin();
95  scIt!= SCs.end(); scIt++){
96 
97  double scEt = scIt->energy() * TMath::Sin(scIt->position().theta());
98  double _eta = scIt->position().eta();
99  if(scEt > scEtMin && _eta > scEtaMin && _eta < scEtaMax) nSC++;
100  }
101  }
102 
103  // Make final desicion on passed events
104  bool accepted = false;
105  if(nJets >= minNumberOfjets || nSC >=minNumberOfSuperClusters) {
106  accepted = true;
107  nSelectedEvents++;
108  }
109 
110  return accepted;
111 }
T getParameter(std::string const &) const
bool filter(edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
int iEvent
Definition: GenABIO.cc:230
EgammaProbeSelector(const edm::ParameterSet &)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
vector< PseudoJet > jets
bool isValid() const
Definition: HandleBase.h:74
T const * product() const
Definition: Handle.h:81
fixed size matrix
HLT enums.
UInt_t nEvents
Definition: hcalCalib.cc:42
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects