CMS 3D CMS Logo

HiggsTo2GammaSkim.cc
Go to the documentation of this file.
1 
2 /* \class HiggsTo2GammaSkim
3  *
4  * Consult header file for description
5  *
6  * author: Kati Lassila-Perini Helsinki Institute of Physics
7  *
8  */
9 
10 // system include files
12 
13 // User include files
15 
16 // Message logger
18 
19 // Photons:
21 
22 // C++
23 #include <iostream>
24 #include <vector>
25 
26 using namespace std;
27 using namespace edm;
28 using namespace reco;
29 
30 // Constructor
32  // Local Debug flag
33  debug = pset.getParameter<bool>("DebugHiggsTo2GammaSkim");
34 
35  // Reconstructed objects
36  thePhotonToken = consumes<reco::PhotonCollection>(pset.getParameter<edm::InputTag>("PhotonCollectionLabel"));
37 
38  // Minimum Pt for photons for skimming
39  photon1MinPt = pset.getParameter<double>("photon1MinimumPt");
40  nPhotonMin = pset.getParameter<int>("nPhotonMinimum");
41 
42  nEvents = 0;
43  nSelectedEvents = 0;
44 }
45 
46 // Destructor
48  edm::LogVerbatim("HiggsTo2GammaSkim") << " Number_events_read " << nEvents << " Number_events_kept "
49  << nSelectedEvents << " Efficiency "
50  << ((double)nSelectedEvents) / ((double)nEvents + 0.01) << std::endl;
51 }
52 
53 // Filter event
55  nEvents++;
56 
58 
59  bool keepEvent = false;
60  int nPhotons = 0;
61 
62  // Look at photons:
63 
64  // Get the photon collection from the event
66 
67  event.getByToken(thePhotonToken, photonHandle);
68 
69  if (photonHandle.isValid()) {
70  const reco::PhotonCollection* phoCollection = photonHandle.product();
71 
72  reco::PhotonCollection::const_iterator photons;
73 
74  // Loop over photon collections and count how many photons there are,
75  // and how many are above the thresholds
76 
77  // Question: do we need to take the reconstructed primary vertex at this point?
78  // Here, I assume that the et is taken with respect to the nominal vertex (0,0,0).
79  for (photons = phoCollection->begin(); photons != phoCollection->end(); ++photons) {
80  float et_p = photons->et();
81  if (et_p > photon1MinPt)
82  nPhotons++;
83  }
84  }
85 
86  // Make decision:
87  if (nPhotons >= nPhotonMin)
88  keepEvent = true;
89 
90  if (keepEvent)
91  nSelectedEvents++;
92 
93  return keepEvent;
94 }
T getParameter(std::string const &) const
~HiggsTo2GammaSkim() override
bool isValid() const
Definition: HandleBase.h:70
#define debug
Definition: HDRShower.cc:19
T const * product() const
Definition: Handle.h:69
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
HiggsTo2GammaSkim(const edm::ParameterSet &)
fixed size matrix
bool filter(edm::Event &, const edm::EventSetup &) override
Get event properties to send to builder to fill seed collection.
HLT enums.
UInt_t nEvents
Definition: hcalCalib.cc:41
Definition: event.py:1