CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PartonJetCorrectionExample Class Reference
Inheritance diagram for PartonJetCorrectionExample:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
 PartonJetCorrectionExample (const edm::ParameterSet &fParameters)
 
virtual ~PartonJetCorrectionExample ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

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.

38  : mInput (fConfig.getParameter <edm::InputTag> ("src")),
39  m_gJ_CorrectorName (fConfig.getParameter <std::string> ("gJetCorrector")),
40  m_qJ_CorrectorName (fConfig.getParameter <std::string> ("qJetCorrector")),
41  m_bJ_CorrectorName (fConfig.getParameter <std::string> ("bJetCorrector")),
42  m_bT_CorrectorName (fConfig.getParameter <std::string> ("bTopCorrector"))
43 {}
virtual PartonJetCorrectionExample::~PartonJetCorrectionExample ( )
inlinevirtual

Definition at line 14 of file PartonJetCorrectionExample.cc.

14 {}

Member Function Documentation

void PartonJetCorrectionExample::analyze ( const edm::Event fEvent,
const edm::EventSetup fSetup 
)
overridevirtual

Implements edm::EDAnalyzer.

Definition at line 45 of file PartonJetCorrectionExample.cc.

References reco::JetCorrector::correction(), mvaPFMET_cff::corrector, gather_cfg::cout, reco::LeafCandidate::eta(), edm::Event::getByLabel(), JetCorrector::getJetCorrector(), metsig::jet, fwrapper::jets, m_bJ_CorrectorName, m_bT_CorrectorName, m_gJ_CorrectorName, m_qJ_CorrectorName, mInput, reco::LeafCandidate::phi(), and reco::LeafCandidate::pt().

45  {
46  // get all correctors
47  const JetCorrector* gJetCorrector = JetCorrector::getJetCorrector (m_gJ_CorrectorName, fSetup);
48  const JetCorrector* qJetCorrector = JetCorrector::getJetCorrector (m_qJ_CorrectorName, fSetup);
49  const JetCorrector* bJetCorrector = JetCorrector::getJetCorrector (m_bJ_CorrectorName, fSetup);
50  const JetCorrector* bTopCorrector = JetCorrector::getJetCorrector (m_bT_CorrectorName, fSetup);
51  const JetCorrector* corrector = 0;
52 
53  // get input jets (supposed to be MC corrected already)
55  fEvent.getByLabel (mInput, jets);
56  // loop over jets
57  for (unsigned ijet = 0; ijet < jets->size(); ++ijet) {
58  const CaloJet& jet = (*jets)[ijet];
59  std::cout << "PartonJetCorrectionExample::analize-> jet #" << ijet;
60  if (ijet%4 == 0) { // assume it is gluon from diJet
61  std::cout << ": use gJ corrections" << std::endl;
62  corrector = gJetCorrector;
63  }
64  else if (ijet%4 == 1) { // assume it is light quark from diJet
65  std::cout << ": use qJ corrections" << std::endl;
66  corrector = qJetCorrector;
67  }
68  else if (ijet%4 == 2) { // assume it is b quark from diJet
69  std::cout << ": use bJ corrections" << std::endl;
70  corrector = bJetCorrector;
71  }
72  else { // assume it is b quark from ttbar
73  std::cout << ": use bT corrections" << std::endl;
74  corrector = bTopCorrector;
75  }
76  // get selected correction for the jet
77  double correction = corrector->correction (jet);
78  // dump it
79  std::cout << " jet pt/eta/phi: " << jet.pt() << '/' << jet.eta() << '/' << jet.phi()
80  << " -> correction factor: " << correction
81  << ", corrected pt: " << jet.pt()*correction
82  << std::endl;
83  }
84 }
Jets made from CaloTowers.
Definition: CaloJet.h:29
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
tuple corrector
Definition: mvaPFMET_cff.py:88
vector< PseudoJet > jets
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
static const JetCorrector * getJetCorrector(const std::string &fName, const edm::EventSetup &fSetup)
retrieve corrector from the event setup. troughs exception if something is missing ...
Definition: JetCorrector.cc:50
tuple cout
Definition: gather_cfg.py:121
virtual double phi() const
momentum azimuthal angle

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().