CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PatBasicAnalyzer.cc
Go to the documentation of this file.
1 #include <map>
2 #include <string>
3 
4 #include "TH1.h"
5 
12 
19 
21 public:
23  explicit PatBasicAnalyzer(const edm::ParameterSet&);
25  ~PatBasicAnalyzer() override;
26 
27 private:
29  void beginJob() override;
31  void analyze(const edm::Event&, const edm::EventSetup&) override;
33  void endJob() override;
34 
35  // simple map to contain all histograms;
36  // histograms are booked in the beginJob()
37  // method
38  std::map<std::string, TH1F*> histContainer_;
39  // plot number of towers per jet
40  TH1F* jetTowers_;
41 
42  // input tokens
49 };
50 
52  : histContainer_(),
53  photonSrcToken_(consumes<edm::View<pat::Photon>>(iConfig.getUntrackedParameter<edm::InputTag>("photonSrc"))),
54  elecSrcToken_(consumes<edm::View<pat::Electron>>(iConfig.getUntrackedParameter<edm::InputTag>("electronSrc"))),
55  muonSrcToken_(consumes<edm::View<pat::Muon>>(iConfig.getUntrackedParameter<edm::InputTag>("muonSrc"))),
56  tauSrcToken_(consumes<edm::View<pat::Tau>>(iConfig.getUntrackedParameter<edm::InputTag>("tauSrc"))),
57  jetSrcToken_(consumes<edm::View<pat::Jet>>(iConfig.getUntrackedParameter<edm::InputTag>("jetSrc"))),
58  metSrcToken_(consumes<edm::View<pat::MET>>(iConfig.getUntrackedParameter<edm::InputTag>("metSrc"))) {}
59 
61 
63  // get electron collection
65  iEvent.getByToken(elecSrcToken_, electrons);
66 
67  // get muon collection
69  iEvent.getByToken(muonSrcToken_, muons);
70 
71  // get tau collection
73  iEvent.getByToken(tauSrcToken_, taus);
74 
75  // get jet collection
77  iEvent.getByToken(jetSrcToken_, jets);
78 
79  // get met collection
81  iEvent.getByToken(metSrcToken_, mets);
82 
83  // get photon collection
85  iEvent.getByToken(photonSrcToken_, photons);
86 
87  // loop over jets
88  size_t nJets = 0;
89  for (edm::View<pat::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); ++jet) {
90  if (jet->pt() > 50) {
91  ++nJets;
92  }
93  // uncomment the following line to fill the
94  // jetTowers_ histogram
95  // jetTowers_->Fill(jet->getCaloConstituents().size());
96  }
97  histContainer_["jets"]->Fill(nJets);
98 
99  // do something similar for the other candidates
100  histContainer_["photons"]->Fill(photons->size());
101  histContainer_["elecs"]->Fill(electrons->size());
102  histContainer_["muons"]->Fill(muons->size());
103  histContainer_["taus"]->Fill(taus->size());
104  histContainer_["met"]->Fill(mets->empty() ? 0 : (*mets)[0].et());
105 }
106 
108  // register to the TFileService
110 
111  // book histograms:
112  // uncomment the following line to book the jetTowers_ histogram
113  //jetTowers_= fs->make<TH1F>("jetTowers", "towers per jet", 90, 0, 90);
114  histContainer_["photons"] = fs->make<TH1F>("photons", "photon multiplicity", 10, 0, 10);
115  histContainer_["elecs"] = fs->make<TH1F>("elecs", "electron multiplicity", 10, 0, 10);
116  histContainer_["muons"] = fs->make<TH1F>("muons", "muon multiplicity", 10, 0, 10);
117  histContainer_["taus"] = fs->make<TH1F>("taus", "tau multiplicity", 10, 0, 10);
118  histContainer_["jets"] = fs->make<TH1F>("jets", "jet multiplicity", 10, 0, 10);
119  histContainer_["met"] = fs->make<TH1F>("met", "missing E_{T}", 20, 0, 100);
120 }
121 
123 
edm::EDGetTokenT< edm::View< pat::Tau > > tauSrcToken_
constexpr char Photon[]
Definition: modules.cc:14
edm::EDGetTokenT< edm::View< pat::Jet > > jetSrcToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void analyze(const edm::Event &, const edm::EventSetup &) override
everything that needs to be done during the event loop
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
std::map< std::string, TH1F * > histContainer_
PatBasicAnalyzer(const edm::ParameterSet &)
default constructor
edm::EDGetTokenT< edm::View< pat::Muon > > muonSrcToken_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< edm::View< pat::MET > > metSrcToken_
vector< PseudoJet > jets
edm::EDGetTokenT< edm::View< pat::Electron > > elecSrcToken_
constexpr char Jet[]
Definition: modules.cc:9
tuple muons
Definition: patZpeak.py:39
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
~PatBasicAnalyzer() override
default destructor
void beginJob() override
everything that needs to be done before the event loop
void endJob() override
everything that needs to be done after the event loop
constexpr char Electron[]
Definition: modules.cc:12
edm::EDGetTokenT< edm::View< pat::Photon > > photonSrcToken_