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