CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QcdHighPtDQM.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author S. Bolognesi, Eric - CERN
5  */
6 
8 
13 
17 
20 #include <vector>
21 
22 #include <string>
23 #include <cmath>
24 
25 using namespace std;
26 using namespace edm;
27 using namespace reco;
28 using namespace math;
29 
30 // Get Jets and MET (no MET plots yet pending converging w/JetMET group)
31 
33  : jetToken_(consumes<CaloJetCollection>(
34  iConfig.getUntrackedParameter<edm::InputTag>("jetTag"))),
35  metToken1_(consumes<CaloMETCollection>(
36  iConfig.getUntrackedParameter<edm::InputTag>("metTag1"))),
37  metToken2_(consumes<CaloMETCollection>(
38  iConfig.getUntrackedParameter<edm::InputTag>("metTag2"))),
39  metToken3_(consumes<CaloMETCollection>(
40  iConfig.getUntrackedParameter<edm::InputTag>("metTag3"))),
41  metToken4_(consumes<CaloMETCollection>(
42  iConfig.getUntrackedParameter<edm::InputTag>("metTag4"))) {
43 }
44 
46 
48  edm::Run const &,
49  edm::EventSetup const &) {
50  iBooker.setCurrentFolder("Physics/QcdHighPt");
51 
52  MEcontainer_["dijet_mass"] = iBooker.book1D(
53  "dijet_mass", "dijet resonance invariant mass, barrel region", 100, 0,
54  1000);
55  MEcontainer_["njets"] =
56  iBooker.book1D("njets", "jet multiplicity", 10, 0, 10);
57  MEcontainer_["etaphi"] = iBooker.book2D("etaphi", "eta/phi distribution", 83,
58  -42, 42, 72, -M_PI, M_PI);
59  MEcontainer_["njets30"] =
60  iBooker.book1D("njets30", "jet multiplicity, pt > 30 GeV", 10, 0, 10);
61 
62  // book histograms for inclusive jet quantities
63  MEcontainer_["inclusive_jet_pt"] = iBooker.book1D(
64  "inclusive_jet_pt", "inclusive jet Pt spectrum", 100, 0, 1000);
65  MEcontainer_["inclusive_jet_pt_barrel"] = iBooker.book1D(
66  "inclusive_jet_pt_barrel", "inclusive jet Pt, eta < 1.3", 100, 0, 1000);
67  MEcontainer_["inclusive_jet_pt_forward"] =
68  iBooker.book1D("inclusive_jet_pt_forward",
69  "inclusive jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
70  MEcontainer_["inclusive_jet_pt_endcap"] =
71  iBooker.book1D("inclusive_jet_pt_endcap",
72  "inclusive jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
73 
74  // book histograms for leading jet quantities
75  MEcontainer_["leading_jet_pt"] =
76  iBooker.book1D("leading_jet_pt", "leading jet Pt", 100, 0, 1000);
77  MEcontainer_["leading_jet_pt_barrel"] = iBooker.book1D(
78  "leading_jet_pt_barrel", "leading jet Pt, eta < 1.3", 100, 0, 1000);
79  MEcontainer_["leading_jet_pt_forward"] =
80  iBooker.book1D("leading_jet_pt_forward",
81  "leading jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
82  MEcontainer_["leading_jet_pt_endcap"] = iBooker.book1D(
83  "leading_jet_pt_endcap", "leading jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
84 
85  // book histograms for met over sum et and met over leading jet pt for various
86  // flavors of MET
87  MEcontainer_["movers_met"] = iBooker.book1D(
88  "movers_met", "MET over Sum ET for basic MET collection", 50, 0, 1);
89  MEcontainer_["moverl_met"] = iBooker.book1D(
90  "moverl_met", "MET over leading jet Pt for basic MET collection", 50, 0,
91  2);
92 
93  MEcontainer_["movers_metho"] = iBooker.book1D(
94  "movers_metho", "MET over Sum ET for MET HO collection", 50, 0, 1);
95  MEcontainer_["moverl_metho"] =
96  iBooker.book1D("moverl_metho",
97  "MET over leading jet Pt for MET HO collection", 50, 0, 2);
98 
99  MEcontainer_["movers_metnohf"] = iBooker.book1D(
100  "movers_metnohf", "MET over Sum ET for MET no HF collection", 50, 0, 1);
101  MEcontainer_["moverl_metnohf"] = iBooker.book1D(
102  "moverl_metnohf", "MET over leading jet Pt for MET no HF collection", 50,
103  0, 2);
104 
105  MEcontainer_["movers_metnohfho"] =
106  iBooker.book1D("movers_metnohfho",
107  "MET over Sum ET for MET no HF HO collection", 50, 0, 1);
108  MEcontainer_["moverl_metnohfho"] = iBooker.book1D(
109  "moverl_metnohfho", "MET over leading jet Pt for MET no HF HO collection",
110  50, 0, 2);
111 
112  // book histograms for EMF fraction for all jets and first 3 jets
113  MEcontainer_["inclusive_jet_EMF"] =
114  iBooker.book1D("inclusive_jet_EMF", "inclusive jet EMF", 50, -1, 1);
115  MEcontainer_["leading_jet_EMF"] =
116  iBooker.book1D("leading_jet_EMF", "leading jet EMF", 50, -1, 1);
117  MEcontainer_["second_jet_EMF"] =
118  iBooker.book1D("second_jet_EMF", "second jet EMF", 50, -1, 1);
119  MEcontainer_["third_jet_EMF"] =
120  iBooker.book1D("third_jet_EMF", "third jet EMF", 50, -1, 1);
121 }
122 
123 // method to calculate MET over Sum ET from a particular MET collection
124 float QcdHighPtDQM::movers(const CaloMETCollection &metcollection) {
125  float metovers = 0;
126  CaloMETCollection::const_iterator met_iter;
127  for (met_iter = metcollection.begin(); met_iter != metcollection.end();
128  ++met_iter) {
129  float mex = met_iter->momentum().x();
130  float mey = met_iter->momentum().y();
131  float met = sqrt(mex * mex + mey * mey);
132  float sumet = met_iter->sumEt();
133  metovers = (met / sumet);
134  }
135  return metovers;
136 }
137 
138 // method to calculate MET over Leading jet PT for a particular MET collection
139 float QcdHighPtDQM::moverl(const CaloMETCollection &metcollection,
140  float &ljpt) {
141  float metoverl = 0;
142  CaloMETCollection::const_iterator met_iter;
143  for (met_iter = metcollection.begin(); met_iter != metcollection.end();
144  ++met_iter) {
145  float mex = met_iter->momentum().x();
146  float mey = met_iter->momentum().y();
147  float met = sqrt(mex * mex + mey * mey);
148  metoverl = (met / ljpt);
149  }
150  return metoverl;
151 }
152 
153 void QcdHighPtDQM::analyze(const Event &iEvent, const EventSetup &iSetup) {
154  // Get Jets
156  iEvent.getByToken(jetToken_, jetHandle);
157  const CaloJetCollection &jets = *jetHandle;
158  CaloJetCollection::const_iterator jet_iter;
159 
160  // Get MET collections
162  iEvent.getByToken(metToken1_, metHandle);
163  const CaloMETCollection &met = *metHandle;
164 
165  edm::Handle<CaloMETCollection> metHOHandle;
166  iEvent.getByToken(metToken2_, metHOHandle);
167  const CaloMETCollection &metHO = *metHOHandle;
168 
169  edm::Handle<CaloMETCollection> metNoHFHandle;
170  iEvent.getByToken(metToken3_, metNoHFHandle);
171  const CaloMETCollection &metNoHF = *metNoHFHandle;
172 
173  edm::Handle<CaloMETCollection> metNoHFHOHandle;
174  iEvent.getByToken(metToken4_, metNoHFHOHandle);
175  const CaloMETCollection &metNoHFHO = *metNoHFHOHandle;
176 
177  // initialize leading jet value and jet multiplicity counter
178  int njets = 0;
179  int njets30 = 0;
180  float leading_jetpt = 0;
181  float leading_jeteta = 0;
182 
183  // initialize variables for picking out leading 2 barrel jets
184  reco::CaloJet leadingbarreljet;
185  reco::CaloJet secondbarreljet;
186  int nbarreljets = 0;
187 
188  // get bins in eta.
189  // Bins correspond to calotower regions.
190 
191  const float etabins[83] = {
192  -5.191, -4.889, -4.716, -4.538, -4.363, -4.191, -4.013, -3.839, -3.664,
193  -3.489, -3.314, -3.139, -2.964, -2.853, -2.650, -2.500, -2.322, -2.172,
194  -2.043, -1.930, -1.830, -1.740, -1.653, -1.566, -1.479, -1.392, -1.305,
195  -1.218, -1.131, -1.044, -.957, -.879, -.783, -.696, -.609, -.522,
196  -.435, -.348, -.261, -.174, -.087, 0, .087, .174, .261,
197  .348, .435, .522, .609, .696, .783, .879, .957, 1.044,
198  1.131, 1.218, 1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830,
199  1.930, 2.043, 2.172, 2.322, 2.500, 2.650, 2.853, 2.964, 3.139,
200  3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.889,
201  5.191};
202 
203  for (jet_iter = jets.begin(); jet_iter != jets.end(); ++jet_iter) {
204  njets++;
205 
206  // get Jet stats
207  float jet_pt = jet_iter->pt();
208  float jet_eta = jet_iter->eta();
209  float jet_phi = jet_iter->phi();
210 
211  // fill jet Pt and jet EMF
212  MEcontainer_["inclusive_jet_pt"]->Fill(jet_pt);
213  MEcontainer_["inclusive_jet_EMF"]->Fill(jet_iter->emEnergyFraction());
214 
215  // pick out up to the first 2 leading barrel jets
216  // for use in calculating dijet mass in barrel region
217  // also fill jet Pt histogram for barrel
218 
219  if (jet_eta <= 1.3) {
220  MEcontainer_["inclusive_jet_pt_barrel"]->Fill(jet_pt);
221  if (nbarreljets == 0) {
222  leadingbarreljet = jets[(njets - 1)];
223  nbarreljets++;
224  } else if (nbarreljets == 1) {
225  secondbarreljet = jets[(njets - 1)];
226  nbarreljets++;
227  }
228 
229  }
230 
231  // fill jet Pt for endcap and forward regions
232  else if (jet_eta <= 3.0 && jet_eta > 1.3) {
233  MEcontainer_["inclusive_jet_pt_endcap"]->Fill(jet_pt);
234  } else if (jet_eta <= 5.0 && jet_eta > 3.0) {
235  MEcontainer_["inclusive_jet_pt_forward"]->Fill(jet_pt);
236  }
237 
238  // count jet multiplicity for jets with Pt > 30
239  if ((jet_pt) > 30) njets30++;
240 
241  // check leading jet quantities
242  if (jet_pt > leading_jetpt) {
243  leading_jetpt = jet_pt;
244  leading_jeteta = jet_eta;
245  }
246 
247  // fill eta-phi plot
248  for (int eit = 0; eit < 81; eit++) {
249  for (int pit = 0; pit < 72; pit++) {
250  float low_eta = etabins[eit];
251  float high_eta = etabins[eit + 1];
252  float low_phi = (-M_PI) + pit * (M_PI / 36);
253  float high_phi = low_phi + (M_PI / 36);
254  if (jet_eta > low_eta && jet_eta < high_eta && jet_phi > low_phi &&
255  jet_phi < high_phi) {
256  MEcontainer_["etaphi"]->Fill((eit - 41), jet_phi);
257  }
258  }
259  }
260  }
261 
262  // after iterating over all jets, fill leading jet quantity histograms
263  // and jet multiplicity histograms
264 
265  MEcontainer_["leading_jet_pt"]->Fill(leading_jetpt);
266 
267  if (leading_jeteta <= 1.3) {
268  MEcontainer_["leading_jet_pt_barrel"]->Fill(leading_jetpt);
269  } else if (leading_jeteta <= 3.0 && leading_jeteta > 1.3) {
270  MEcontainer_["leading_jet_pt_endcap"]->Fill(leading_jetpt);
271  } else if (leading_jeteta <= 5.0 && leading_jeteta > 3.0) {
272  MEcontainer_["leading_jet_pt_forward"]->Fill(leading_jetpt);
273  }
274 
275  MEcontainer_["njets"]->Fill(njets);
276 
277  MEcontainer_["njets30"]->Fill(njets30);
278 
279  // fill MET over Sum ET and Leading jet PT for all MET flavors
280  MEcontainer_["movers_met"]->Fill(movers(met));
281  MEcontainer_["moverl_met"]->Fill(movers(met), leading_jetpt);
282  MEcontainer_["movers_metho"]->Fill(movers(metHO));
283  MEcontainer_["moverl_metho"]->Fill(movers(metHO), leading_jetpt);
284  MEcontainer_["movers_metnohf"]->Fill(movers(metNoHF));
285  MEcontainer_["moverl_metnohf"]->Fill(movers(metNoHF), leading_jetpt);
286  MEcontainer_["movers_metnohfho"]->Fill(movers(metNoHFHO));
287  MEcontainer_["moverl_metnohfho"]->Fill(movers(metNoHFHO), leading_jetpt);
288 
289  // fetch first 3 jet EMF
290 
291  if (jets.size() >= 1) {
292  MEcontainer_["leading_jet_EMF"]->Fill(jets[0].emEnergyFraction());
293  if (jets.size() >= 2) {
294  MEcontainer_["second_jet_EMF"]->Fill(jets[1].emEnergyFraction());
295  if (jets.size() >= 3) {
296  MEcontainer_["third_jet_EMF"]->Fill(jets[2].emEnergyFraction());
297  }
298  }
299  }
300 
301  // if 2 nontrivial barrel jets, reconstruct dijet mass
302 
303  if (nbarreljets == 2) {
304  if (leadingbarreljet.energy() > 0 && secondbarreljet.energy() > 0) {
306  leadingbarreljet.p4() + secondbarreljet.p4();
307  float dijet_mass = DiJet.mass();
308  MEcontainer_["dijet_mass"]->Fill(dijet_mass);
309  }
310  }
311 }
312 
313 // Local Variables:
314 // show-trailing-whitespace: t
315 // truncate-lines: t
316 // End:
edm::EDGetTokenT< reco::CaloMETCollection > metToken1_
Definition: QcdHighPtDQM.h:34
float moverl(const reco::CaloMETCollection &metcollection, float &ljpt)
void analyze(const edm::Event &, const edm::EventSetup &)
Jets made from CaloTowers.
Definition: CaloJet.h:29
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
edm::EDGetTokenT< reco::CaloMETCollection > metToken2_
Definition: QcdHighPtDQM.h:35
virtual double pt() const
transverse momentum
std::map< std::string, MonitorElement * > MEcontainer_
Definition: QcdHighPtDQM.h:40
Collection of Calo MET.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
virtual double energy() const
energy
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:48
vector< PseudoJet > jets
edm::EDGetTokenT< reco::CaloMETCollection > metToken4_
Definition: QcdHighPtDQM.h:37
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
edm::EDGetTokenT< reco::CaloMETCollection > metToken3_
Definition: QcdHighPtDQM.h:36
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: QcdHighPtDQM.cc:47
#define M_PI
float movers(const reco::CaloMETCollection &metcollection)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
QcdHighPtDQM(const edm::ParameterSet &)
Definition: QcdHighPtDQM.cc:32
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
edm::EDGetTokenT< reco::CaloJetCollection > jetToken_
Definition: QcdHighPtDQM.h:33
virtual ~QcdHighPtDQM()
Definition: QcdHighPtDQM.cc:45
Definition: Run.h:41
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects