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
JetPtAnalyzer Class Reference

#include <JetPtAnalyzer.h>

Inheritance diagram for JetPtAnalyzer:
JetAnalyzerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &, const reco::CaloJetCollection &caloJets)
 Get the analysis. More...
 
void beginJob (DQMStore *dbe)
 Inizialize parameters for histo binning. More...
 
void endJob ()
 Finish up a job. More...
 
 JetPtAnalyzer (const edm::ParameterSet &)
 Constructor. More...
 
void setSource (std::string source)
 
virtual ~JetPtAnalyzer ()
 Destructor. More...
 
- Public Member Functions inherited from JetAnalyzerBase
void analyze (const edm::Event &, const edm::EventSetup &, reco::CaloJet &caloJet)
 Get the analysis of the muon properties. More...
 
 JetAnalyzerBase ()
 Constructor. More...
 
virtual ~JetAnalyzerBase ()
 Destructor. More...
 

Private Attributes

std::string _source
 
int etaBin
 
double etaMax
 
double etaMin
 
reco::helper::JetIDHelperjetID
 
MonitorElementjetME
 
std::string jetname
 
MonitorElementmConstituents
 
MonitorElementmEFrac
 
MonitorElementmEmEnergyInEB
 
MonitorElementmEmEnergyInEE
 
MonitorElementmEmEnergyInHF
 
MonitorElementmEta
 
MonitorElementmfHPD
 
MonitorElementmfRBX
 
MonitorElementmHadEnergyInHB
 
MonitorElementmHadEnergyInHE
 
MonitorElementmHadEnergyInHF
 
MonitorElementmHadEnergyInHO
 
MonitorElementmHFrac
 
MonitorElementmMaxEInEmTowers
 
MonitorElementmMaxEInHadTowers
 
MonitorElementmN90Hits
 
MonitorElementmNJets
 
MonitorElementmPhi
 
MonitorElementmresEMF
 
edm::ParameterSet parameters
 
int phiBin
 
double phiMax
 
double phiMin
 
int ptBin
 
double ptMax
 
double ptMin
 
edm::InputTag theCaloJetCollectionLabel
 

Detailed Description

DQM monitoring source for Calo Jets

Author
J.Weng ETH Zuerich

Definition at line 35 of file JetPtAnalyzer.h.

Constructor & Destructor Documentation

JetPtAnalyzer::JetPtAnalyzer ( const edm::ParameterSet pSet)

Constructor.

Definition at line 19 of file JetPtAnalyzer.cc.

References Parameters::parameters.

19  {
20 
21  parameters = pSet;
22 
23  // Note: one could also add cut on energies from short and long fibers
24  // or fraction of energy from the hottest HPD readout
25 
26  //-------------------
27 }
edm::ParameterSet parameters
Definition: JetPtAnalyzer.h:62
JetPtAnalyzer::~JetPtAnalyzer ( )
virtual

Destructor.

Definition at line 30 of file JetPtAnalyzer.cc.

30 { }

Member Function Documentation

void JetPtAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const reco::CaloJetCollection caloJets 
)

Get the analysis.

Definition at line 95 of file JetPtAnalyzer.cc.

References metsig::jet, bTagSequences_cff::jetID, and LogTrace.

96  {
97  // int numofjets=0;
98  LogTrace(jetname)<<"[JetPtAnalyzer] Analyze Calo Jet";
99 
100  for (reco::CaloJetCollection::const_iterator jet = caloJets.begin(); jet!=caloJets.end(); ++jet){
101  //-----------------------------
102  jetME->Fill(1);
103 
104  jetID->calculate(iEvent, *jet);
105 
106  if (mEta) mEta->Fill (jet->pt(),jet->eta());
107  if (mPhi) mPhi->Fill (jet->pt(),jet->phi());
108 
109  // if (mNJets) mNJets->Fill (jet->pt(),numofjets);
110  if (mConstituents) mConstituents->Fill (jet->pt(),jet->nConstituents());
111  if (mHFrac) mHFrac->Fill (jet->pt(),jet->energyFractionHadronic());
112  if (mEFrac) mEFrac->Fill (jet->pt(),jet->emEnergyFraction());
113 
114 
115  if (mMaxEInEmTowers) mMaxEInEmTowers->Fill (jet->pt(),jet->maxEInEmTowers());
116  if (mMaxEInHadTowers) mMaxEInHadTowers->Fill (jet->pt(),jet->maxEInHadTowers());
117 
118  if (mHadEnergyInHO) mHadEnergyInHO->Fill (jet->pt(),jet->hadEnergyInHO());
119  if (mHadEnergyInHB) mHadEnergyInHB->Fill (jet->pt(),jet->hadEnergyInHB());
120  if (mHadEnergyInHF) mHadEnergyInHF->Fill (jet->pt(),jet->hadEnergyInHF());
121  if (mHadEnergyInHE) mHadEnergyInHE->Fill (jet->pt(),jet->hadEnergyInHE());
122  if (mEmEnergyInEB) mEmEnergyInEB->Fill (jet->pt(),jet->emEnergyInEB());
123  if (mEmEnergyInEE) mEmEnergyInEE->Fill (jet->pt(),jet->emEnergyInEE());
124  if (mEmEnergyInHF) mEmEnergyInHF->Fill (jet->pt(),jet->emEnergyInHF());
125 
126  if (mN90Hits) mN90Hits->Fill (jet->pt(),jetID->n90Hits());
127  if (mresEMF) mresEMF->Fill(jet->pt(),jetID->restrictedEMF());
128  if (mfHPD) mfHPD->Fill(jet->pt(),jetID->fHPD());
129  if (mfRBX) mfRBX->Fill(jet->pt(),jetID->fRBX());
130 
131  }
132 }
MonitorElement * mEmEnergyInEB
MonitorElement * mHadEnergyInHF
MonitorElement * mN90Hits
MonitorElement * mConstituents
Definition: JetPtAnalyzer.h:90
MonitorElement * mMaxEInEmTowers
Definition: JetPtAnalyzer.h:97
MonitorElement * mHadEnergyInHO
Definition: JetPtAnalyzer.h:99
double fHPD() const
Definition: JetIDHelper.h:34
double restrictedEMF() const
Definition: JetIDHelper.h:51
MonitorElement * mHadEnergyInHB
MonitorElement * mfRBX
MonitorElement * mresEMF
MonitorElement * jetME
Definition: JetPtAnalyzer.h:86
std::string jetname
Definition: JetPtAnalyzer.h:64
void Fill(long long x)
MonitorElement * mEFrac
Definition: JetPtAnalyzer.h:92
MonitorElement * mHadEnergyInHE
#define LogTrace(id)
MonitorElement * mEmEnergyInEE
MonitorElement * mEmEnergyInHF
MonitorElement * mEta
Definition: JetPtAnalyzer.h:87
MonitorElement * mMaxEInHadTowers
Definition: JetPtAnalyzer.h:98
MonitorElement * mHFrac
Definition: JetPtAnalyzer.h:91
MonitorElement * mfHPD
MonitorElement * mPhi
Definition: JetPtAnalyzer.h:88
reco::helper::JetIDHelper * jetID
Definition: JetPtAnalyzer.h:70
void calculate(const edm::Event &event, const reco::CaloJet &jet, const int iDbg=0)
Definition: JetIDHelper.cc:86
double fRBX() const
Definition: JetIDHelper.h:35
void JetPtAnalyzer::beginJob ( DQMStore dbe)
virtual

Inizialize parameters for histo binning.

Implements JetAnalyzerBase.

Definition at line 34 of file JetPtAnalyzer.cc.

References DQMStore::book1D(), DQMStore::book2D(), jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, bTagSequences_cff::jetID, LogTrace, Parameters::parameters, jptDQMConfig_cff::phiMax, jptDQMConfig_cff::phiMin, jptDQMConfig_cff::ptMax, PtMinSelector_cfg::ptMin, MonitorElement::setBinLabel(), and DQMStore::setCurrentFolder().

34  {
35 
36  jetname = "jetPtAnalyzer";
37 
38  LogTrace(jetname)<<"[JetPtAnalyzer] Parameters initialization";
39  dbe->setCurrentFolder("JetMET/Jet/"+_source);
40 
41  jetME = dbe->book1D("jetReco", "jetReco", 3, 1, 4);
42  jetME->setBinLabel(1,"CaloJets",1);
43 
44  //initialize JetID
46 
47  // monitoring of eta parameter
48  etaBin = parameters.getParameter<int>("etaBin");
49  etaMin = parameters.getParameter<double>("etaMin");
50  etaMax = parameters.getParameter<double>("etaMax");
51 
52  // monitoring of phi paramater
53  phiBin = parameters.getParameter<int>("phiBin");
54  phiMin = parameters.getParameter<double>("phiMin");
55  phiMax = parameters.getParameter<double>("phiMax");
56 
57  // monitoring of the transverse momentum
58  ptBin = parameters.getParameter<int>("ptBin");
59  ptMin = parameters.getParameter<double>("ptMin");
60  ptMax = parameters.getParameter<double>("ptMax");
61 
62  mEta = dbe->book2D("EtaVsPt", "EtaVsPt",ptBin, ptMin, ptMax, etaBin, etaMin, etaMax);
63  mPhi = dbe->book2D("PhiVsPt", "PhiVsPt",ptBin, ptMin, ptMax, phiBin, phiMin, phiMax);
64  mConstituents = dbe->book2D("ConstituentsVsPt", "# of ConstituentsVsPt",ptBin, ptMin, ptMax, 100, 0, 100);
65  // mNJets = dbe->book2D("NJetsVsPt", "number of jetsVsPt",ptBin, ptMin, ptMax, 100, 0, 100);
66 
67  // CaloJet specific
68  mHFrac = dbe->book2D("HFracVsPt", "HFracVsPt",ptBin, ptMin, ptMax, 120, -0.1, 1.1);
69  mEFrac = dbe->book2D("EFracVsPt", "EFracVsPt", ptBin, ptMin, ptMax,120, -0.1, 1.1);
70  mMaxEInEmTowers = dbe->book2D("MaxEInEmTowersVsPt", "MaxEInEmTowersVsPt", ptBin, ptMin, ptMax,100, 0, 100);
71  mMaxEInHadTowers = dbe->book2D("MaxEInHadTowersVsPt", "MaxEInHadTowersVsPt",ptBin, ptMin, ptMax, 100, 0, 100);
72  mHadEnergyInHO = dbe->book2D("HadEnergyInHOVsPt", "HadEnergyInHOVsPt",ptBin, ptMin, ptMax, 100, 0, 10);
73  mHadEnergyInHB = dbe->book2D("HadEnergyInHBVsPt", "HadEnergyInHBVsPt",ptBin, ptMin, ptMax, 100, 0, 50);
74  mHadEnergyInHF = dbe->book2D("HadEnergyInHFVsPt", "HadEnergyInHFVsPt",ptBin, ptMin, ptMax ,100, 0, 50);
75  mHadEnergyInHE = dbe->book2D("HadEnergyInHEVsPt", "HadEnergyInHEVsPt",ptBin, ptMin, ptMax, 100, 0, 100);
76  mEmEnergyInEB = dbe->book2D("EmEnergyInEBVsPt", "EmEnergyInEBVsPt",ptBin, ptMin, ptMax, 100, 0, 50);
77  mEmEnergyInEE = dbe->book2D("EmEnergyInEEVsPt", "EmEnergyInEEVsPt",ptBin, ptMin, ptMax, 100, 0, 50);
78  mEmEnergyInHF = dbe->book2D("EmEnergyInHFVsPt", "EmEnergyInHFVsPt",ptBin, ptMin, ptMax, 120, -20, 100);
79  //jetID variables
80  mN90Hits = dbe->book2D("N90HitsVsPt", "N90HitsVsPt",ptBin, ptMin, ptMax, 50, 0, 50);
81  mresEMF = dbe->book2D("resEMFVsPt", "resEMFVsPt",ptBin, ptMin, ptMax,50, 0., 1.);
82  mfHPD = dbe->book2D("fHPDVsPt", "fHPDVsPt",ptBin, ptMin, ptMax,50, 0., 1.);
83  mfRBX = dbe->book2D("fRBXVsPt", "fRBXVsPt",ptBin, ptMin, ptMax,50, 0., 1.);
84 
85 
86 }
MonitorElement * mEmEnergyInEB
T getParameter(std::string const &) const
MonitorElement * mHadEnergyInHF
MonitorElement * mN90Hits
MonitorElement * mConstituents
Definition: JetPtAnalyzer.h:90
MonitorElement * mMaxEInEmTowers
Definition: JetPtAnalyzer.h:97
edm::ParameterSet parameters
Definition: JetPtAnalyzer.h:62
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
MonitorElement * mHadEnergyInHO
Definition: JetPtAnalyzer.h:99
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * mHadEnergyInHB
MonitorElement * mfRBX
MonitorElement * mresEMF
MonitorElement * jetME
Definition: JetPtAnalyzer.h:86
std::string jetname
Definition: JetPtAnalyzer.h:64
MonitorElement * mEFrac
Definition: JetPtAnalyzer.h:92
MonitorElement * mHadEnergyInHE
#define LogTrace(id)
MonitorElement * mEmEnergyInEE
MonitorElement * mEmEnergyInHF
MonitorElement * mEta
Definition: JetPtAnalyzer.h:87
MonitorElement * mMaxEInHadTowers
Definition: JetPtAnalyzer.h:98
std::string _source
Definition: JetPtAnalyzer.h:65
MonitorElement * mHFrac
Definition: JetPtAnalyzer.h:91
MonitorElement * mfHPD
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1000
MonitorElement * mPhi
Definition: JetPtAnalyzer.h:88
reco::helper::JetIDHelper * jetID
Definition: JetPtAnalyzer.h:70
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
void JetPtAnalyzer::endJob ( void  )

Finish up a job.

Definition at line 87 of file JetPtAnalyzer.cc.

References bTagSequences_cff::jetID.

87  {
88 
89  delete jetID;
90 
91 }
reco::helper::JetIDHelper * jetID
Definition: JetPtAnalyzer.h:70
void JetPtAnalyzer::setSource ( std::string  source)
inline

Definition at line 55 of file JetPtAnalyzer.h.

References _source, and source.

55  {
56  _source = source;
57  }
std::string _source
Definition: JetPtAnalyzer.h:65
static std::string const source
Definition: EdmProvDump.cc:43

Member Data Documentation

std::string JetPtAnalyzer::_source
private

Definition at line 65 of file JetPtAnalyzer.h.

Referenced by setSource().

int JetPtAnalyzer::etaBin
private

Definition at line 73 of file JetPtAnalyzer.h.

double JetPtAnalyzer::etaMax
private

Definition at line 75 of file JetPtAnalyzer.h.

double JetPtAnalyzer::etaMin
private

Definition at line 74 of file JetPtAnalyzer.h.

reco::helper::JetIDHelper* JetPtAnalyzer::jetID
private

Definition at line 70 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::jetME
private

Definition at line 86 of file JetPtAnalyzer.h.

std::string JetPtAnalyzer::jetname
private

Definition at line 64 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mConstituents
private

Definition at line 90 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mEFrac
private

Definition at line 92 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mEmEnergyInEB
private

Definition at line 103 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mEmEnergyInEE
private

Definition at line 104 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mEmEnergyInHF
private

Definition at line 105 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mEta
private

Definition at line 87 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mfHPD
private

Definition at line 108 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mfRBX
private

Definition at line 109 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mHadEnergyInHB
private

Definition at line 100 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mHadEnergyInHE
private

Definition at line 102 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mHadEnergyInHF
private

Definition at line 101 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mHadEnergyInHO
private

Definition at line 99 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mHFrac
private

Definition at line 91 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mMaxEInEmTowers
private

Definition at line 97 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mMaxEInHadTowers
private

Definition at line 98 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mN90Hits
private

Definition at line 106 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mNJets
private

Definition at line 94 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mPhi
private

Definition at line 88 of file JetPtAnalyzer.h.

MonitorElement* JetPtAnalyzer::mresEMF
private

Definition at line 107 of file JetPtAnalyzer.h.

edm::ParameterSet JetPtAnalyzer::parameters
private
int JetPtAnalyzer::phiBin
private

Definition at line 77 of file JetPtAnalyzer.h.

double JetPtAnalyzer::phiMax
private

Definition at line 79 of file JetPtAnalyzer.h.

double JetPtAnalyzer::phiMin
private

Definition at line 78 of file JetPtAnalyzer.h.

int JetPtAnalyzer::ptBin
private

Definition at line 81 of file JetPtAnalyzer.h.

double JetPtAnalyzer::ptMax
private

Definition at line 83 of file JetPtAnalyzer.h.

double JetPtAnalyzer::ptMin
private

Definition at line 82 of file JetPtAnalyzer.h.

edm::InputTag JetPtAnalyzer::theCaloJetCollectionLabel
private

Definition at line 67 of file JetPtAnalyzer.h.