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
FlavorJetCorrectionExample Class Reference
Inheritance diagram for FlavorJetCorrectionExample:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
 FlavorJetCorrectionExample (const edm::ParameterSet &fParameters)
 
virtual ~FlavorJetCorrectionExample ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

std::string mBCorrectorName
 
std::string mCCorrectorName
 
edm::InputTag mInput
 
std::string mUDSCorrectorName
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- 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::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 11 of file FlavorJetCorrectionExample.cc.

Constructor & Destructor Documentation

FlavorJetCorrectionExample::FlavorJetCorrectionExample ( const edm::ParameterSet fParameters)
explicit

Definition at line 36 of file FlavorJetCorrectionExample.cc.

37  : mInput (fConfig.getParameter <edm::InputTag> ("src")),
38  mUDSCorrectorName (fConfig.getParameter <std::string> ("UDSQuarksCorrector")),
39  mCCorrectorName (fConfig.getParameter <std::string> ("CQuarkCorrector")),
40  mBCorrectorName (fConfig.getParameter <std::string> ("BQuarkCorrector"))
41 {}
virtual FlavorJetCorrectionExample::~FlavorJetCorrectionExample ( )
inlinevirtual

Definition at line 14 of file FlavorJetCorrectionExample.cc.

14 {}

Member Function Documentation

void FlavorJetCorrectionExample::analyze ( const edm::Event fEvent,
const edm::EventSetup fSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 43 of file FlavorJetCorrectionExample.cc.

References JetCorrector::correction(), gather_cfg::cout, reco::LeafCandidate::eta(), edm::Event::getByLabel(), JetCorrector::getJetCorrector(), metsig::jet, fwrapper::jets, mBCorrectorName, mCCorrectorName, mInput, mUDSCorrectorName, reco::LeafCandidate::phi(), and reco::LeafCandidate::pt().

43  {
44  // get all correctors
45  const JetCorrector* udsJetCorrector = JetCorrector::getJetCorrector (mUDSCorrectorName, fSetup);
46  const JetCorrector* cQuarkJetCorrector = JetCorrector::getJetCorrector (mCCorrectorName, fSetup);
47  const JetCorrector* bQuarkJetCorrector = JetCorrector::getJetCorrector (mBCorrectorName, fSetup);
48  const JetCorrector* corrector = 0;
49 
50  // get input jets (supposed to be MC corrected already)
52  fEvent.getByLabel (mInput, jets);
53  // loop over jets
54  for (unsigned ijet = 0; ijet < jets->size(); ++ijet) {
55  const CaloJet& jet = (*jets)[ijet];
56  std::cout << "FlavorJetCorrectionExample::analize-> jet #" << ijet;
57  if (ijet%3 == 0) { // assume it is light quark
58  std::cout << ": use USD quark corrections" << std::endl;
59  corrector = udsJetCorrector;
60  }
61  else if (ijet%3 == 1) { // assume it is c quark
62  std::cout << ": use c quark corrections" << std::endl;
63  corrector = cQuarkJetCorrector;
64  }
65  else { // assume it is b quark
66  std::cout << ": use b quark corrections" << std::endl;
67  corrector = bQuarkJetCorrector;
68  }
69  // get selected correction for the jet
70  double correction = corrector->correction (jet);
71  // dump it
72  std::cout << " jet pt/eta/phi: " << jet.pt() << '/' << jet.eta() << '/' << jet.phi()
73  << " -> correction factor: " << correction
74  << ", corrected pt: " << jet.pt()*correction
75  << std::endl;
76  }
77 }
Jets made from CaloTowers.
Definition: CaloJet.h:30
virtual double correction(const LorentzVector &fJet) const =0
get correction using Jet information only
virtual double eta() const
momentum pseudorapidity
vector< PseudoJet > jets
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual double pt() const
transverse momentum
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:51
tuple cout
Definition: gather_cfg.py:121
virtual double phi() const
momentum azimuthal angle

Member Data Documentation

std::string FlavorJetCorrectionExample::mBCorrectorName
private

Definition at line 20 of file FlavorJetCorrectionExample.cc.

Referenced by analyze().

std::string FlavorJetCorrectionExample::mCCorrectorName
private

Definition at line 19 of file FlavorJetCorrectionExample.cc.

Referenced by analyze().

edm::InputTag FlavorJetCorrectionExample::mInput
private

Definition at line 17 of file FlavorJetCorrectionExample.cc.

Referenced by analyze().

std::string FlavorJetCorrectionExample::mUDSCorrectorName
private

Definition at line 18 of file FlavorJetCorrectionExample.cc.

Referenced by analyze().