CMS 3D CMS Logo

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