CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
17 public:
18  explicit PatZjetsJetAnalyzer(const edm::ParameterSet&);
20 
21 private:
22 
23  virtual void beginJob() override ;
24  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
25  virtual void endJob() override ;
26 
27  // simple map to contain all histograms;
28  // histograms are booked in the beginJob()
29  // method
30  std::map<std::string,TH1F*> histContainer_;
31 
32  // input tags
34 };
35 
37  histContainer_(),
38  srcToken_(consumes<edm::View<pat::Jet> >(iConfig.getUntrackedParameter<edm::InputTag>("src")))
39 {
40 }
41 
43 {
44 }
45 
46 void
48 {
49  // get electron collection
51  iEvent.getByToken(srcToken_,jets);
52 
53  // loop jets
54  for(edm::View<pat::Jet>::const_iterator ijet=jets->begin(); ijet!=jets->end(); ++ijet){
55  // fill simple histograms
56  pat::Jet jet = ijet->correctedJet("had", "uds");
57  histContainer_["pt" ]->Fill( jet.pt () );
58  histContainer_["eta" ]->Fill( jet.eta() );
59  histContainer_["phi" ]->Fill( jet.phi() );
60  histContainer_["emf" ]->Fill( jet.emEnergyFraction() );
61  for(unsigned int i=0; i<jet.getCaloConstituents().size(); ++i){
62  histContainer_["dEta"]->Fill( jet.getCaloConstituent(i)->eta()-jet.eta() );
63  }
64  }
65 }
66 
67 void
69 {
70  // register to the TFileService
72 
73  // book histograms:
74  histContainer_["pt" ]=fs->make<TH1F>("pt" , "pt" , 150, 0., 150.);
75  histContainer_["eta" ]=fs->make<TH1F>("eta" , "eta" , 50, 0., 5.);
76  histContainer_["phi" ]=fs->make<TH1F>("phi" , "phi" , 60, 3.14, 3.14);
77  histContainer_["emf" ]=fs->make<TH1F>("emf" , "emf" , 40, 0., 1.);
78  histContainer_["dEta"]=fs->make<TH1F>("dEta" , "dEta" , 40, 0., 1.);
79 }
80 
81 void
83 {
84 }
85 
int i
Definition: DBlmapReader.cc:9
PatZjetsJetAnalyzer(const edm::ParameterSet &)
virtual void endJob() override
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
virtual void beginJob() override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual double phi() const final
momentum azimuthal angle
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:230
vector< PseudoJet > jets
float emEnergyFraction() const
returns the jet electromagnetic energy fraction
Definition: Jet.h:296
edm::EDGetTokenT< edm::View< pat::Jet > > srcToken_
std::map< std::string, TH1F * > histContainer_
Analysis-level calorimeter jet class.
Definition: Jet.h:78
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
virtual double eta() const final
momentum pseudorapidity
CaloTowerPtr getCaloConstituent(unsigned fIndex) const
convert generic constituent to specific type
edm::Service< TFileService > fs
virtual double pt() const final
transverse momentum