CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/TopQuarkAnalysis/Examples/plugins/TopJetAnalyzer.cc

Go to the documentation of this file.
00001 #include "DataFormats/PatCandidates/interface/Jet.h"
00002 #include "TopQuarkAnalysis/Examples/plugins/TopJetAnalyzer.h"
00003 
00004 
00005 TopJetAnalyzer::TopJetAnalyzer(const edm::ParameterSet& cfg):
00006   input_  (cfg.getParameter<edm::InputTag>("input"  )),
00007   verbose_(cfg.getParameter<bool>         ("verbose"))
00008 {
00009   edm::Service<TFileService> fs;
00010   
00011   mult_ = fs->make<TH1F>("mult", "multiplicity (jets)", 30,  0 ,   30);
00012   en_   = fs->make<TH1F>("en"  , "energy (jets)",       60,  0., 300.);
00013   pt_   = fs->make<TH1F>("pt"  , "pt (jets)",           60,  0., 300.);
00014   eta_  = fs->make<TH1F>("eta" , "eta (jets)",          30, -3.,   3.);
00015   phi_  = fs->make<TH1F>("phi" , "phi (jets)",          40, -4.,   4.);
00016 }
00017 
00018 TopJetAnalyzer::~TopJetAnalyzer()
00019 {
00020 }
00021 
00022 void
00023 TopJetAnalyzer::analyze(const edm::Event& evt, const edm::EventSetup& setup)
00024 {
00025   edm::Handle<std::vector<pat::Jet> > jets;
00026   evt.getByLabel(input_, jets); 
00027 
00028   // fill histograms
00029   
00030   mult_->Fill( jets->size() );
00031   for(std::vector<pat::Jet>::const_iterator jet=jets->begin(); jet!=jets->end(); ++jet){
00032     pt_ ->Fill( jet->pt()     );
00033     en_ ->Fill( jet->energy() );
00034     eta_->Fill( jet->eta()    );
00035     phi_->Fill( jet->phi()    );
00036   }
00037 
00038   // produce printout if desired
00039 
00040   if( jets->size()<1 || !verbose_ )
00041     return;
00042 
00043   int lineWidth = 75;
00044   if( jets->begin()->isCaloJet() )
00045     lineWidth = 100;
00046   else if( jets->begin()->isPFJet() )
00047     lineWidth = 120;
00048 
00049   std::cout << std::setfill('=') << std::setw(lineWidth) << "\n" << std::setfill(' ');
00050   std::cout << std::setw( 5) << "jet :"
00051             << std::setw(11) << "pt :"
00052             << std::setw( 9) << "eta :"
00053             << std::setw( 9) << "phi :"
00054             << std::setw(11) << "TCHE :"
00055             << std::setw(11) << "TCHP :"
00056             << std::setw( 9) << "SSVHE :"
00057             << std::setw( 9) << "SSVHP :";
00058   if( jets->begin()->isCaloJet() ) {
00059     std::cout << std::setw( 8) << "emf :"
00060               << std::setw(10) << "n90Hits :"
00061               << std::setw( 7) << "fHPD";
00062   }
00063   if( jets->begin()->isPFJet() ) {
00064     std::cout << std::setw(9) << "chf : "
00065               << std::setw(8) << "nhf : "
00066               << std::setw(8) << "cef : "
00067               << std::setw(8) << "nef : "
00068               << std::setw(6) << "nCh : "
00069               << std::setw(6) << "nConst";
00070   }
00071   std::cout << std::endl
00072             << std::setfill('-') << std::setw(lineWidth) << "\n" << std::setfill(' ');
00073   unsigned i=0;
00074   for(std::vector<pat::Jet>::const_iterator jet=jets->begin(); jet!=jets->end(); ++jet){
00075     std::cout << std::setw(3) << i << " : " << std::setprecision(3) << std::fixed
00076               << std::setw(8) << jet->pt() << " : "
00077               << std::setw(6) << jet->eta() << " : "
00078               << std::setw(6) << jet->phi() << " : "
00079               << std::setw(8) << jet->bDiscriminator("trackCountingHighEffBJetTags") << " : "
00080               << std::setw(8) << jet->bDiscriminator("trackCountingHighPurBJetTags") << " : "
00081               << std::setw(6) << jet->bDiscriminator("simpleSecondaryVertexHighEffBJetTags") << " : "
00082               << std::setw(6) << jet->bDiscriminator("simpleSecondaryVertexHighPurBJetTags") << " : ";
00083     if( jet->isCaloJet() ) {
00084       std::cout << std::setw(5) << jet->emEnergyFraction() << " : "
00085                 << std::setw(7) << jet->jetID().n90Hits << " : "
00086                 << std::setw(6) << jet->jetID().fHPD;
00087     }
00088     if( jet->isPFJet() ) {
00089       std::cout << std::setw(5) << jet->chargedHadronEnergyFraction() << " : "
00090                 << std::setw(5) << jet->neutralHadronEnergyFraction() << " : "
00091                 << std::setw(5) << jet->chargedEmEnergyFraction() << " : "
00092                 << std::setw(5) << jet->neutralEmEnergyFraction() << " : "
00093                 << std::setw(3) << jet->chargedMultiplicity() << " : "
00094                 << std::setw(6) << jet->nConstituents();
00095     }
00096     std::cout << std::endl;
00097     i++;
00098   }
00099   std::cout << std::setfill('=') << std::setw(lineWidth) << "\n" << std::setfill(' ');
00100 }
00101 
00102 void TopJetAnalyzer::beginJob()
00103 {
00104 }
00105 
00106 void TopJetAnalyzer::endJob()
00107 {
00108 }
00109 
00110