CMS 3D CMS Logo

TopJetAnalyzer.cc
Go to the documentation of this file.
2 
4  : inputToken_(consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("input"))),
5  verbose_(cfg.getParameter<bool>("verbose")) {
7 
8  mult_ = fs->make<TH1F>("mult", "multiplicity (jets)", 30, 0, 30);
9  en_ = fs->make<TH1F>("en", "energy (jets)", 60, 0., 300.);
10  pt_ = fs->make<TH1F>("pt", "pt (jets)", 60, 0., 300.);
11  eta_ = fs->make<TH1F>("eta", "eta (jets)", 30, -3., 3.);
12  phi_ = fs->make<TH1F>("phi", "phi (jets)", 40, -4., 4.);
13 }
14 
16 
19  evt.getByToken(inputToken_, jets);
20 
21  // fill histograms
22 
23  mult_->Fill(jets->size());
24  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); ++jet) {
25  pt_->Fill(jet->pt());
26  en_->Fill(jet->energy());
27  eta_->Fill(jet->eta());
28  phi_->Fill(jet->phi());
29  }
30 
31  // produce printout if desired
32 
33  if (jets->empty() || !verbose_)
34  return;
35 
36  int lineWidth = 75;
37  if (jets->begin()->isCaloJet())
38  lineWidth = 100;
39  else if (jets->begin()->isPFJet())
40  lineWidth = 120;
41 
42  std::cout << std::setfill('=') << std::setw(lineWidth) << "\n" << std::setfill(' ');
43  std::cout << std::setw(5) << "jet :" << std::setw(11) << "pt :" << std::setw(9) << "eta :" << std::setw(9)
44  << "phi :" << std::setw(11) << "TCHE :" << std::setw(11) << "TCHP :" << std::setw(9)
45  << "SSVHE :" << std::setw(9) << "SSVHP :";
46  if (jets->begin()->isCaloJet()) {
47  std::cout << std::setw(8) << "emf :" << std::setw(10) << "n90Hits :" << std::setw(7) << "fHPD";
48  }
49  if (jets->begin()->isPFJet()) {
50  std::cout << std::setw(9) << "chf : " << std::setw(8) << "nhf : " << std::setw(8) << "cef : " << std::setw(8)
51  << "nef : " << std::setw(6) << "nCh : " << std::setw(6) << "nConst";
52  }
53  std::cout << std::endl << std::setfill('-') << std::setw(lineWidth) << "\n" << std::setfill(' ');
54  unsigned i = 0;
55  for (std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); ++jet) {
56  std::cout << std::setw(3) << i << " : " << std::setprecision(3) << std::fixed << std::setw(8) << jet->pt() << " : "
57  << std::setw(6) << jet->eta() << " : " << std::setw(6) << jet->phi() << " : " << std::setw(8)
58  << jet->bDiscriminator("trackCountingHighEffBJetTags") << " : " << std::setw(8)
59  << jet->bDiscriminator("trackCountingHighPurBJetTags") << " : " << std::setw(6)
60  << jet->bDiscriminator("simpleSecondaryVertexHighEffBJetTags") << " : " << std::setw(6)
61  << jet->bDiscriminator("simpleSecondaryVertexHighPurBJetTags") << " : ";
62  if (jet->isCaloJet()) {
63  std::cout << std::setw(5) << jet->emEnergyFraction() << " : " << std::setw(7) << jet->jetID().n90Hits << " : "
64  << std::setw(6) << jet->jetID().fHPD;
65  }
66  if (jet->isPFJet()) {
67  std::cout << std::setw(5) << jet->chargedHadronEnergyFraction() << " : " << std::setw(5)
68  << jet->neutralHadronEnergyFraction() << " : " << std::setw(5) << jet->chargedEmEnergyFraction()
69  << " : " << std::setw(5) << jet->neutralEmEnergyFraction() << " : " << std::setw(3)
70  << jet->chargedMultiplicity() << " : " << std::setw(6) << jet->nConstituents();
71  }
72  std::cout << std::endl;
73  i++;
74  }
75  std::cout << std::setfill('=') << std::setw(lineWidth) << "\n" << std::setfill(' ');
76 }
77 
79 
~TopJetAnalyzer() override
void beginJob() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
Definition: HeavyIon.h:7
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: Jet.py:1
void endJob() override
TopJetAnalyzer(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< pat::Jet > > inputToken_
HLT enums.