CMS 3D CMS Logo

PartonJetCorrectionExample Class Reference

Inheritance diagram for PartonJetCorrectionExample:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 PartonJetCorrectionExample (const edm::ParameterSet &fParameters)
virtual ~PartonJetCorrectionExample ()

Private Attributes

std::string m_bJ_CorrectorName
std::string m_bT_CorrectorName
std::string m_gJ_CorrectorName
std::string m_qJ_CorrectorName
edm::InputTag mInput


Detailed Description

Definition at line 11 of file PartonJetCorrectionExample.cc.


Constructor & Destructor Documentation

PartonJetCorrectionExample::PartonJetCorrectionExample ( const edm::ParameterSet fParameters  )  [explicit]

Definition at line 37 of file PartonJetCorrectionExample.cc.

00038   : mInput (fConfig.getParameter <edm::InputTag> ("src")),
00039     m_gJ_CorrectorName (fConfig.getParameter <std::string> ("gJetCorrector")),
00040     m_qJ_CorrectorName (fConfig.getParameter <std::string> ("qJetCorrector")),
00041     m_bJ_CorrectorName (fConfig.getParameter <std::string> ("bJetCorrector")),
00042     m_bT_CorrectorName (fConfig.getParameter <std::string> ("bTopCorrector"))
00043 {}

virtual PartonJetCorrectionExample::~PartonJetCorrectionExample (  )  [inline, virtual]

Definition at line 14 of file PartonJetCorrectionExample.cc.

00014 {}


Member Function Documentation

void PartonJetCorrectionExample::analyze ( const edm::Event fEvent,
const edm::EventSetup fSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 45 of file PartonJetCorrectionExample.cc.

References JetCorrector::correction(), GenMuonPlsPt100GeV_cfg::cout, lat::endl(), reco::Particle::eta(), edm::Event::getByLabel(), JetCorrector::getJetCorrector(), metsig::jet, pfTauBenchmarkGeneric_cfi::jets, m_bJ_CorrectorName, m_bT_CorrectorName, m_gJ_CorrectorName, m_qJ_CorrectorName, mInput, reco::Particle::phi(), and reco::Particle::pt().

00045                                                                                             {
00046   // get all correctors
00047   const JetCorrector* gJetCorrector = JetCorrector::getJetCorrector (m_gJ_CorrectorName, fSetup);
00048   const JetCorrector* qJetCorrector = JetCorrector::getJetCorrector (m_qJ_CorrectorName, fSetup);
00049   const JetCorrector* bJetCorrector = JetCorrector::getJetCorrector (m_bJ_CorrectorName, fSetup);
00050   const JetCorrector* bTopCorrector = JetCorrector::getJetCorrector (m_bT_CorrectorName, fSetup);
00051   const JetCorrector* corrector = 0;
00052   
00053   // get input jets (supposed to be MC corrected already)
00054   edm::Handle<CaloJetCollection> jets;                    
00055   fEvent.getByLabel (mInput, jets);
00056   // loop over jets                      
00057   for (unsigned ijet = 0; ijet < jets->size(); ++ijet) {
00058     const CaloJet& jet = (*jets)[ijet];
00059     std::cout << "PartonJetCorrectionExample::analize-> jet #" << ijet;
00060     if (ijet%4 == 0) { // assume it is gluon from diJet
00061       std::cout << ": use gJ corrections" << std::endl;
00062       corrector = gJetCorrector;
00063     }
00064     else if (ijet%4 == 1) { // assume it is light quark from diJet
00065       std::cout << ": use qJ corrections" << std::endl;
00066       corrector = qJetCorrector;
00067     }
00068     else if (ijet%4 == 2) { // assume it is b quark from diJet
00069       std::cout << ": use bJ corrections" << std::endl;
00070       corrector = bJetCorrector;
00071     }
00072     else { // assume it is b quark from ttbar
00073       std::cout << ": use bT corrections" << std::endl;
00074       corrector = bTopCorrector;
00075     }
00076     // get selected correction for the jet
00077     double correction = corrector->correction (jet, fEvent, fSetup);
00078     // dump it
00079     std::cout << "  jet pt/eta/phi: " << jet.pt() << '/' <<  jet.eta() << '/' << jet.phi() 
00080               << " -> correction factor: " << correction 
00081               << ", corrected pt: " << jet.pt()*correction
00082               << std::endl;
00083   }
00084 }


Member Data Documentation

std::string PartonJetCorrectionExample::m_bJ_CorrectorName [private]

Definition at line 20 of file PartonJetCorrectionExample.cc.

Referenced by analyze().

std::string PartonJetCorrectionExample::m_bT_CorrectorName [private]

Definition at line 21 of file PartonJetCorrectionExample.cc.

Referenced by analyze().

std::string PartonJetCorrectionExample::m_gJ_CorrectorName [private]

Definition at line 18 of file PartonJetCorrectionExample.cc.

Referenced by analyze().

std::string PartonJetCorrectionExample::m_qJ_CorrectorName [private]

Definition at line 19 of file PartonJetCorrectionExample.cc.

Referenced by analyze().

edm::InputTag PartonJetCorrectionExample::mInput [private]

Definition at line 17 of file PartonJetCorrectionExample.cc.

Referenced by analyze().


The documentation for this class was generated from the following file:
Generated on Tue Jun 9 18:29:27 2009 for CMSSW by  doxygen 1.5.4