CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetPtAnalyzer.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author J.Weng ETH Zuerich
5  */
6 
9 
12 
14 
15 #include <string>
16 using namespace edm;
17 
18 // ***********************************************************
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 }
28 
29 // ***********************************************************
31 
32 
33 // ***********************************************************
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
45  jetID = new reco::helper::JetIDHelper(parameters.getParameter<ParameterSet>("JetIDParams"));
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 }
88 
89  delete jetID;
90 
91 }
92 
93 
94 // ***********************************************************
96  const reco::CaloJetCollection& caloJets) {
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 }
dictionary parameters
Definition: Parameters.py:2
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
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)
virtual ~JetPtAnalyzer()
Destructor.
int iEvent
Definition: GenABIO.cc:243
void analyze(const edm::Event &, const edm::EventSetup &, const reco::CaloJetCollection &caloJets)
Get the analysis.
#define LogTrace(id)
void beginJob(DQMStore *dbe)
Inizialize parameters for histo binning.
JetPtAnalyzer(const edm::ParameterSet &)
Constructor.
void endJob()
Finish up a job.
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:845
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects