CMS 3D CMS Logo

AlCaGammaJetSelector.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Calibration/HcalAlCaRecoProducers/AlCaGammaJetSelector
4 // Class: AlCaGammaJetSelector
5 //
13 //
14 // Original Author: Andrius Juodagalvis
15 // Created: Fri, 15 Aug 2015 12:09:55 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <iostream>
22 
23 // user include files
29 
33 //#include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
34 
35 //
36 // class declaration
37 //
38 
39 namespace AlCaGammaJet {
40  struct Counters {
42  mutable std::atomic<unsigned int> nProcessed_, nSelected_;
43  };
44 } // namespace AlCaGammaJet
45 
46 class AlCaGammaJetSelector : public edm::stream::EDFilter<edm::GlobalCache<AlCaGammaJet::Counters> > {
47 public:
49  ~AlCaGammaJetSelector() override;
50 
51  static std::unique_ptr<AlCaGammaJet::Counters> initializeGlobalCache(edm::ParameterSet const&) {
52  return std::make_unique<AlCaGammaJet::Counters>();
53  }
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
56  bool filter(edm::Event&, const edm::EventSetup&) override;
57  void endStream() override;
58  static void globalEndJob(const AlCaGammaJet::Counters* counters);
59 
60 private:
61  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
62  void endRun(edm::Run const&, edm::EventSetup const&) override {}
64  void endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override {}
66 
67  // ----------member data ---------------------------
68 
69  unsigned int nProcessed_, nSelected_;
70 
72  double minPtJet_, minPtPhoton_;
75 };
76 
77 //
78 // constants, enums and typedefs
79 //
80 
81 //
82 // static data member definitions
83 //
84 
85 //
86 // constructors and destructor
87 //
89  nProcessed_ = 0;
90  nSelected_ = 0;
91 
92  // get input
93  labelPhoton_ = iConfig.getParameter<edm::InputTag>("PhoInput");
94  labelPFJet_ = iConfig.getParameter<edm::InputTag>("PFjetInput");
95  minPtJet_ = iConfig.getParameter<double>("MinPtJet");
96  minPtPhoton_ = iConfig.getParameter<double>("MinPtPhoton");
97 
98  // Register consumption
99  tok_Photon_ = consumes<reco::PhotonCollection>(labelPhoton_);
100  tok_PFJet_ = consumes<reco::PFJetCollection>(labelPFJet_);
101 }
102 
104 
105 //
106 // member functions
107 //
108 
109 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
111  //The following says we do not know what parameters are allowed so do no validation
112  // Please change this to state exactly what you do use, even if it is no parameters
114  desc.add<edm::InputTag>("PhoInput", edm::InputTag("gedPhotons"));
115  desc.add<edm::InputTag>("PFjetInput", edm::InputTag("ak4PFJetsCHS"));
116  desc.add<double>("MinPtJet", 10.0);
117  desc.add<double>("MinPtPhoton", 10.0);
118  descriptions.addDefault(desc);
119 }
120 
121 // ------------ method called on each new Event ------------
123  nProcessed_++;
124 
125  // Access the collections from iEvent
127  iEvent.getByToken(tok_Photon_, phoHandle);
128  if (!phoHandle.isValid()) {
129  edm::LogWarning("AlCaGammaJet") << "AlCaGammaJetProducer: Error! can't get the product " << labelPhoton_;
130  return false; // do not filter
131  }
132  const reco::PhotonCollection photons = *(phoHandle.product());
133 
135  iEvent.getByToken(tok_PFJet_, pfjetHandle);
136  if (!pfjetHandle.isValid()) {
137  edm::LogWarning("AlCaGammaJet") << "AlCaGammaJetProducer: Error! can't get product " << labelPFJet_;
138  return false; // do not filter
139  }
140  const reco::PFJetCollection pfjets = *(pfjetHandle.product());
141 
142  // Check the conditions for a good event
143  if (!(select(photons, pfjets)))
144  return false;
145 
146  //std::cout << "good event\n";
147  nSelected_++;
148  return true;
149 }
150 
152  globalCache()->nProcessed_ += nProcessed_;
153  globalCache()->nSelected_ += nSelected_;
154 }
155 
156 // ------------ method called once each job just after ending the event loop ------------
157 
159  edm::LogWarning("AlCaGammaJet") << "Finds " << count->nSelected_ << " good events out of " << count->nProcessed_;
160 }
161 
163  // Check the requirement for minimum pT
164  if (photons.empty())
165  return false;
166  bool ok(false);
167  for (reco::PFJetCollection::const_iterator itr = jets.begin(); itr != jets.end(); ++itr) {
168  if (itr->pt() >= minPtJet_) {
169  ok = true;
170  break;
171  }
172  }
173  if (!ok)
174  return ok;
175  for (reco::PhotonCollection::const_iterator itr = photons.begin(); itr != photons.end(); ++itr) {
176  if (itr->pt() >= minPtPhoton_)
177  return ok;
178  }
179  return false;
180 }
181 
182 //define this as a plug-in
T getParameter(std::string const &) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
static void globalEndJob(const AlCaGammaJet::Counters *counters)
void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
std::atomic< unsigned int > nSelected_
AlCaGammaJetSelector(const edm::ParameterSet &, const AlCaGammaJet::Counters *count)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
bool filter(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< reco::PFJetCollection > tok_PFJet_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void endRun(edm::Run const &, edm::EventSetup const &) override
T const * product() const
Definition: Handle.h:69
bool select(const reco::PhotonCollection &, const reco::PFJetCollection &)
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
static std::unique_ptr< AlCaGammaJet::Counters > initializeGlobalCache(edm::ParameterSet const &)
void beginRun(edm::Run const &, edm::EventSetup const &) override
select
when omitted electron plots will be filled w/o cut on electronId electronId = cms.PSet( src = cms.InputTag("mvaTrigV0"), cutValue = cms.double(0.5) ), when omitted electron plots will be filled w/o additional pre- selection of the electron candidates
std::vector< PFJet > PFJetCollection
collection of PFJet objects
std::atomic< unsigned int > nProcessed_
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
edm::EDGetTokenT< reco::PhotonCollection > tok_Photon_
Definition: Run.h:45