CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PatZjetsJetAnalyzer.cc
Go to the documentation of this file.
1 #include <map>
2 #include <string>
3 
4 #include "TH1.h"
5 
12 
14 
16 public:
17  explicit PatZjetsJetAnalyzer(const edm::ParameterSet&);
18  ~PatZjetsJetAnalyzer() override;
19 
20 private:
21  void beginJob() override;
22  void analyze(const edm::Event&, const edm::EventSetup&) override;
23  void endJob() override;
24 
25  // simple map to contain all histograms;
26  // histograms are booked in the beginJob()
27  // method
28  std::map<std::string, TH1F*> histContainer_;
29 
30  // input tags
32 };
33 
35  : histContainer_(),
36  srcToken_(consumes<edm::View<pat::Jet> >(iConfig.getUntrackedParameter<edm::InputTag>("src"))) {}
37 
39 
41  // get electron collection
43  iEvent.getByToken(srcToken_, jets);
44 
45  // loop jets
46  for (edm::View<pat::Jet>::const_iterator ijet = jets->begin(); ijet != jets->end(); ++ijet) {
47  // fill simple histograms
48  pat::Jet jet = ijet->correctedJet("had", "uds");
49  histContainer_["pt"]->Fill(jet.pt());
50  histContainer_["eta"]->Fill(jet.eta());
51  histContainer_["phi"]->Fill(jet.phi());
52  histContainer_["emf"]->Fill(jet.emEnergyFraction());
53  for (unsigned int i = 0; i < jet.getCaloConstituents().size(); ++i) {
54  histContainer_["dEta"]->Fill(jet.getCaloConstituent(i)->eta() - jet.eta());
55  }
56  }
57 }
58 
60  // register to the TFileService
62 
63  // book histograms:
64  histContainer_["pt"] = fs->make<TH1F>("pt", "pt", 150, 0., 150.);
65  histContainer_["eta"] = fs->make<TH1F>("eta", "eta", 50, 0., 5.);
66  histContainer_["phi"] = fs->make<TH1F>("phi", "phi", 60, 3.14, 3.14);
67  histContainer_["emf"] = fs->make<TH1F>("emf", "emf", 40, 0., 1.);
68  histContainer_["dEta"] = fs->make<TH1F>("dEta", "dEta", 40, 0., 1.);
69 }
70 
72 
PatZjetsJetAnalyzer(const edm::ParameterSet &)
double pt() const final
transverse momentum
void analyze(const edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void beginJob() override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
std::vector< CaloTowerPtr > const & getCaloConstituents() const
int iEvent
Definition: GenABIO.cc:224
vector< PseudoJet > jets
float emEnergyFraction() const
returns the jet electromagnetic energy fraction
Definition: Jet.h:326
edm::EDGetTokenT< edm::View< pat::Jet > > srcToken_
Analysis-level calorimeter jet class.
Definition: Jet.h:77
constexpr char Jet[]
Definition: modules.cc:9
std::map< std::string, TH1F * > histContainer_
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
CaloTowerPtr getCaloConstituent(unsigned fIndex) const
convert generic constituent to specific type
double phi() const final
momentum azimuthal angle
double eta() const final
momentum pseudorapidity