#include <TopJetAnalyzer.h>
Definition at line 13 of file TopJetAnalyzer.h.
Definition at line 5 of file TopJetAnalyzer.cc.
References en_, eta_, TFileDirectory::make(), mult_, phi_, and pt_.
11 mult_ = fs->
make<TH1F>(
"mult",
"multiplicity (jets)", 30, 0 , 30);
12 en_ = fs->
make<TH1F>(
"en" ,
"energy (jets)", 60, 0., 300.);
13 pt_ = fs->make<TH1F>(
"pt" ,
"pt (jets)", 60, 0., 300.);
14 eta_ = fs->make<TH1F>(
"eta" ,
"eta (jets)", 30, -3., 3.);
15 phi_ = fs->make<TH1F>(
"phi" ,
"phi (jets)", 40, -4., 4.);
T getParameter(std::string const &) const
T * make() const
make new ROOT object
TopJetAnalyzer::~TopJetAnalyzer |
( |
| ) |
|
Implements edm::EDAnalyzer.
Definition at line 23 of file TopJetAnalyzer.cc.
References gather_cfg::cout, en_, eta_, edm::Event::getByLabel(), i, input_, metsig::jet, fwrapper::jets, mult_, phi_, pt_, and verbose_.
30 mult_->Fill( jets->size() );
31 for(std::vector<pat::Jet>::const_iterator
jet=jets->begin();
jet!=jets->end(); ++
jet){
33 en_ ->Fill(
jet->energy() );
44 if( jets->begin()->isCaloJet() )
46 else if( jets->begin()->isPFJet() )
49 std::cout << std::setfill(
'=') << std::setw(lineWidth) <<
"\n" << std::setfill(
' ');
51 << std::setw(11) <<
"pt :"
52 << std::setw( 9) <<
"eta :"
53 << std::setw( 9) <<
"phi :"
54 << std::setw(11) <<
"TCHE :"
55 << std::setw(11) <<
"TCHP :"
56 << std::setw( 9) <<
"SSVHE :"
57 << std::setw( 9) <<
"SSVHP :";
58 if( jets->begin()->isCaloJet() ) {
60 << std::setw(10) <<
"n90Hits :"
61 << std::setw( 7) <<
"fHPD";
63 if( jets->begin()->isPFJet() ) {
65 << std::setw(8) <<
"nhf : "
66 << std::setw(8) <<
"cef : "
67 << std::setw(8) <<
"nef : "
68 << std::setw(6) <<
"nCh : "
69 << std::setw(6) <<
"nConst";
72 << std::setfill(
'-') << std::setw(lineWidth) <<
"\n" << std::setfill(
' ');
74 for(std::vector<pat::Jet>::const_iterator
jet=jets->begin();
jet!=jets->end(); ++
jet){
75 std::cout << std::setw(3) << i <<
" : " << std::setprecision(3) << std::fixed
76 << std::setw(8) <<
jet->pt() <<
" : "
77 << std::setw(6) <<
jet->eta() <<
" : "
78 << std::setw(6) <<
jet->phi() <<
" : "
79 << std::setw(8) <<
jet->bDiscriminator(
"trackCountingHighEffBJetTags") <<
" : "
80 << std::setw(8) <<
jet->bDiscriminator(
"trackCountingHighPurBJetTags") <<
" : "
81 << std::setw(6) <<
jet->bDiscriminator(
"simpleSecondaryVertexHighEffBJetTags") <<
" : "
82 << std::setw(6) <<
jet->bDiscriminator(
"simpleSecondaryVertexHighPurBJetTags") <<
" : ";
83 if(
jet->isCaloJet() ) {
84 std::cout << std::setw(5) <<
jet->emEnergyFraction() <<
" : "
85 << std::setw(7) <<
jet->jetID().n90Hits <<
" : "
86 << std::setw(6) <<
jet->jetID().fHPD;
88 if(
jet->isPFJet() ) {
89 std::cout << std::setw(5) <<
jet->chargedHadronEnergyFraction() <<
" : "
90 << std::setw(5) <<
jet->neutralHadronEnergyFraction() <<
" : "
91 << std::setw(5) <<
jet->chargedEmEnergyFraction() <<
" : "
92 << std::setw(5) <<
jet->neutralEmEnergyFraction() <<
" : "
93 << std::setw(3) <<
jet->chargedMultiplicity() <<
" : "
94 << std::setw(6) <<
jet->nConstituents();
99 std::cout << std::setfill(
'=') << std::setw(lineWidth) <<
"\n" << std::setfill(
' ');
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
void TopJetAnalyzer::beginJob |
( |
void |
| ) |
|
|
privatevirtual |
void TopJetAnalyzer::endJob |
( |
void |
| ) |
|
|
privatevirtual |
TH1F* TopJetAnalyzer::en_ |
|
private |
TH1F* TopJetAnalyzer::eta_ |
|
private |
TH1F* TopJetAnalyzer::mult_ |
|
private |
TH1F* TopJetAnalyzer::phi_ |
|
private |
TH1F* TopJetAnalyzer::pt_ |
|
private |
bool TopJetAnalyzer::verbose_ |
|
private |