CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
22 public:
24  explicit PatBasicAnalyzer(const edm::ParameterSet&);
27 
28 private:
30  virtual void beginJob() override ;
32  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
34  virtual void endJob() override ;
35 
36  // simple map to contain all histograms;
37  // histograms are booked in the beginJob()
38  // method
39  std::map<std::string,TH1F*> histContainer_;
40  // plot number of towers per jet
41  TH1F* jetTowers_;
42 
43  // input tokens
50 };
51 
53  histContainer_(),
54  photonSrcToken_(consumes<edm::View<pat::Photon> >(iConfig.getUntrackedParameter<edm::InputTag>("photonSrc"))),
55  elecSrcToken_(consumes<edm::View<pat::Electron> >(iConfig.getUntrackedParameter<edm::InputTag>("electronSrc"))),
56  muonSrcToken_(consumes<edm::View<pat::Muon> >(iConfig.getUntrackedParameter<edm::InputTag>("muonSrc"))),
57  tauSrcToken_(consumes<edm::View<pat::Tau> >(iConfig.getUntrackedParameter<edm::InputTag>("tauSrc"))),
58  jetSrcToken_(consumes<edm::View<pat::Jet> >(iConfig.getUntrackedParameter<edm::InputTag>("jetSrc"))),
59  metSrcToken_(consumes<edm::View<pat::MET> >(iConfig.getUntrackedParameter<edm::InputTag>("metSrc")))
60 {
61 }
62 
64 {
65 }
66 
67 void
69 {
70  // get electron collection
72  iEvent.getByToken(elecSrcToken_,electrons);
73 
74  // get muon collection
76  iEvent.getByToken(muonSrcToken_,muons);
77 
78  // get tau collection
80  iEvent.getByToken(tauSrcToken_,taus);
81 
82  // get jet collection
84  iEvent.getByToken(jetSrcToken_,jets);
85 
86  // get met collection
88  iEvent.getByToken(metSrcToken_,mets);
89 
90  // get photon collection
92  iEvent.getByToken(photonSrcToken_,photons);
93 
94  // loop over jets
95  size_t nJets=0;
96  for(edm::View<pat::Jet>::const_iterator jet=jets->begin(); jet!=jets->end(); ++jet){
97  if(jet->pt()>50){
98  ++nJets;
99  }
100  // uncomment the following line to fill the
101  // jetTowers_ histogram
102  // jetTowers_->Fill(jet->getCaloConstituents().size());
103  }
104  histContainer_["jets"]->Fill(nJets);
105 
106  // do something similar for the other candidates
107  histContainer_["photons"]->Fill(photons->size() );
108  histContainer_["elecs" ]->Fill(electrons->size());
109  histContainer_["muons"]->Fill(muons->size() );
110  histContainer_["taus" ]->Fill(taus->size() );
111  histContainer_["met" ]->Fill(mets->empty() ? 0 : (*mets)[0].et());
112 }
113 
114 void
116 {
117  // register to the TFileService
119 
120  // book histograms:
121  // uncomment the following line to book the jetTowers_ histogram
122  //jetTowers_= fs->make<TH1F>("jetTowers", "towers per jet", 90, 0, 90);
123  histContainer_["photons"]=fs->make<TH1F>("photons", "photon multiplicity", 10, 0, 10);
124  histContainer_["elecs" ]=fs->make<TH1F>("elecs", "electron multiplicity", 10, 0, 10);
125  histContainer_["muons" ]=fs->make<TH1F>("muons", "muon multiplicity", 10, 0, 10);
126  histContainer_["taus" ]=fs->make<TH1F>("taus", "tau multiplicity", 10, 0, 10);
127  histContainer_["jets" ]=fs->make<TH1F>("jets", "jet multiplicity", 10, 0, 10);
128  histContainer_["met" ]=fs->make<TH1F>("met", "missing E_{T}", 20, 0, 100);
129 }
130 
131 void
133 {
134 }
135 
edm::EDGetTokenT< edm::View< pat::Tau > > tauSrcToken_
edm::EDGetTokenT< edm::View< pat::Jet > > jetSrcToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual 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
PatBasicAnalyzer(const edm::ParameterSet &)
default constructor
edm::EDGetTokenT< edm::View< pat::Muon > > muonSrcToken_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::View< pat::MET > > metSrcToken_
vector< PseudoJet > jets
~PatBasicAnalyzer()
default destructor
std::map< std::string, TH1F * > histContainer_
edm::EDGetTokenT< edm::View< pat::Electron > > elecSrcToken_
tuple muons
Definition: patZpeak.py:38
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
virtual void beginJob() override
everything that needs to be done before the event loop
virtual void endJob() override
everything that needs to be done after the event loop
edm::EDGetTokenT< edm::View< pat::Photon > > photonSrcToken_