CMS 3D CMS Logo

Lepjets_Event.cc
Go to the documentation of this file.
1 //
2 //
3 // File: src/Lepjets_Event.h
4 // Purpose: Represent a simple `event' consisting of leptons and jets.
5 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
6 //
7 // CMSSW File : src/Lepjets_Event.cc
8 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
9 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
10 //
11 
12 
36 #include <algorithm>
37 #include <functional>
38 #include <cmath>
39 #include <cassert>
40 
41 
42 using std::stable_sort;
43 using std::remove_if;
44 using std::abs;
45 
46 
47 namespace hitfit {
48 
49 
51 //
52 // Purpose: Constructor.
53 //
54 // Inputs:
55 // runnum - The run number.
56 // evnum - The event number.
57 //
58  : _zvertex (0),
59  _isMC(false),
60  _runnum (runnum),
61  _evnum (evnum),
62  _dlb (-1),
63  _dnn (-1)
64 {
65 }
66 
67 
68 const int& Lepjets_Event::runnum () const
69 //
70 // Purpose: Return the run number.
71 //
72 // Returns:
73 // The run number.
74 //
75 {
76  return _runnum;
77 }
78 
79 
81 //
82 // Purpose: Return the run number.
83 //
84 // Returns:
85 // The run number.
86 //
87 {
88  return _runnum;
89 }
90 
91 
92 const int& Lepjets_Event::evnum () const
93 //
94 // Purpose: Return the event number.
95 //
96 // Returns:
97 // The event number.
98 //
99 {
100  return _evnum;
101 }
102 
103 
105 //
106 // Purpose: Return the event number.
107 //
108 // Returns:
109 // The event number.
110 //
111 {
112  return _evnum;
113 }
114 
115 
117 //
118 // Purpose: Return the length of the lepton list.
119 //
120 // Returns:
121 // The length of the lepton list.
122 //
123 {
124  return _leps.size ();
125 }
126 
127 
129 //
130 // Purpose: Return the length of the jet list.
131 //
132 // Returns:
133 // The length of the jet list.
134 //
135 {
136  return _jets.size ();
137 }
138 
139 
141 //
142 // Purpose: Return the Ith lepton.
143 //
144 // Inputs:
145 // i - The lepton index (counting from 0).
146 //
147 // Returns:
148 // The Ith lepton.
149 //
150 {
151  assert (i < _leps.size());
152  return _leps[i];
153 }
154 
155 
157 //
158 // Purpose: Return the Ith jet.
159 //
160 // Inputs:
161 // i - The jet index (counting from 0).
162 //
163 // Returns:
164 // The Ith jet.
165 //
166 {
167  assert (i < _jets.size());
168  return _jets[i];
169 }
170 
171 
173 //
174 // Purpose: Return the Ith lepton.
175 //
176 // Inputs:
177 // i - The lepton index (counting from 0).
178 //
179 // Returns:
180 // The Ith lepton.
181 //
182 {
183  assert (i < _leps.size());
184  return _leps[i];
185 }
186 
187 
189 //
190 // Purpose: Return the Ith jet.
191 //
192 // Inputs:
193 // i - The jet index (counting from 0).
194 //
195 // Returns:
196 // The Ith jet.
197 //
198 {
199  assert (i < _jets.size());
200  return _jets[i];
201 }
202 
203 
205 //
206 // Purpose: Return the missing Et.
207 //
208 // Returns:
209 // The missing Et.
210 //
211 {
212  return _met;
213 }
214 
215 
217 //
218 // Purpose: Return the missing Et.
219 //
220 // Returns:
221 // The missing Et.
222 //
223 {
224  return _met;
225 }
226 
227 
229 //
230 // Purpose: Return the kt resolution.
231 //
232 // Returns:
233 // The kt resolution.
234 //
235 {
236  return _kt_res;
237 }
238 
239 
241 //
242 // Purpose: Return the kt resolution.
243 //
244 // Returns:
245 // The kt resolution.
246 //
247 {
248  return _kt_res;
249 }
250 
251 
252 double Lepjets_Event::zvertex () const
253 //
254 // Purpose: Return the z-vertex.
255 //
256 // Returns:
257 // The z-vertex.
258 //
259 {
260  return _zvertex;
261 }
262 
263 
265 //
266 // Purpose: Return the z-vertex.
267 //
268 // Returns:
269 // The z-vertex.
270 //
271 {
272  return _zvertex;
273 }
274 
275 
276 bool Lepjets_Event::isMC () const
277 //
278 // Purpose: Return the isMC flag.
279 //
280 // Returns:
281 // The isMC flag.
282 //
283 {
284  return _isMC;
285 }
286 
287 
289 //
290 // Purpose: set isMC flag.
291 //
292 // Returns:
293 // nothing
294 //
295 {
296  _isMC = isMC;
297 }
298 
299 double Lepjets_Event::dlb () const
300 //
301 // Purpose: Return the LB discriminant.
302 //
303 // Returns:
304 // The LB discriminant.
305 //
306 {
307  return _dlb;
308 }
309 
310 
312 //
313 // Purpose: Return the LB discriminant.
314 //
315 // Returns:
316 // The LB discriminant.
317 //
318 {
319  return _dlb;
320 }
321 
322 
323 double Lepjets_Event::dnn () const
324 //
325 // Purpose: Return the NN discriminant.
326 //
327 // Returns:
328 // The NN discriminant.
329 //
330 {
331  return _dnn;
332 }
333 
334 
336 //
337 // Purpose: Return the NN discriminant.
338 //
339 // Returns:
340 // The NN discriminant.
341 //
342 {
343  return _dnn;
344 }
345 
346 
348 //
349 // Purpose: Sum all objects with type code TYPE.
350 //
351 // Inputs:
352 // type - The type code to match.
353 //
354 // Returns:
355 // The sum of all objects with type code TYPE.
356 //
357 {
358  Fourvec out;
360  if (_leps[i].type() == type)
361  out += _leps[i].p();
363  if (_jets[i].type() == type)
364  out += _jets[i].p();
365  return out;
366 }
367 
368 
370 //
371 // Purpose: Calculate kt --- sum of all objects plus missing Et.
372 //
373 // Returns:
374 // The event kt.
375 {
376  Fourvec v = _met;
378  v += _leps[i].p();
380  v += _jets[i].p();
381  return v;
382 }
383 
384 
386 //
387 // Purpose: Add a lepton to the event.
388 //
389 // Inputs:
390 // lep - The lepton to add.
391 //
392 {
393  _leps.push_back (lep);
394 }
395 
396 
398 //
399 // Purpose: Add a jet to the event.
400 //
401 // Inputs:
402 // jet - The jet to add.
403 //
404 {
405  _jets.push_back (jet);
406 }
407 
408 
409 void Lepjets_Event::smear (CLHEP::HepRandomEngine& engine, bool smear_dir /*= false*/)
410 //
411 // Purpose: Smear the objects in the event according to their resolutions.
412 //
413 // Inputs:
414 // engine - The underlying RNG.
415 // smear_dir - If false, smear the momentum only.
416 //
417 {
418  Fourvec before, after;
419  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
420  before += _leps[i].p();
421  _leps[i].smear (engine, smear_dir);
422  after += _leps[i].p();
423  }
424  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
425  before += _jets[i].p();
426  _jets[i].smear (engine, smear_dir);
427  after += _jets[i].p();
428  }
429 
430  Fourvec kt = _met + before;
433  _met = kt - after;
434 }
435 
436 
438 //
439 // Purpose: Sort the objects in the event in order of descending pt.
440 //
441 {
442  std::stable_sort (_leps.begin(), _leps.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
443  std::stable_sort (_jets.begin(), _jets.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
444 }
445 
446 
447 std::vector<int> Lepjets_Event::jet_types() const
448 //
449 // Purpose: Return the jet types of the event
450 //
451 {
452  std::vector<int> ret;
454  ijet != _jets.size() ;
455  ijet++) {
456  ret.push_back(jet(ijet).type());
457  }
458  return ret;
459 }
460 
461 
462 bool Lepjets_Event::set_jet_types(const std::vector<int>& _jet_types)
463 //
464 // Purpose: Set the jet types of the event
465 // Return false if it fails, trus if it succeeds
466 //
467 {
468  if (_jets.size() != _jet_types.size()) {
469  return false;
470  }
471  bool saw_hadw1 = false;
473  int t = _jet_types[i];
474  if (t == hadw1_label) {
475  if (saw_hadw1)
476  t = hadw2_label;
477  saw_hadw1 = true;
478  }
479  jet (i).type() = t;
480  }
481  return true;
482 }
483 
484 
485 namespace {
486 
487 
488 struct Lepjets_Event_Cutter
489 //
490 // Purpose: Helper for cutting on objects.
491 //
492 {
493  Lepjets_Event_Cutter (double pt_cut, double eta_cut)
494  : _pt_cut (pt_cut), _eta_cut (eta_cut)
495  {}
496  bool operator() (const Lepjets_Event_Lep& l) const;
497  double _pt_cut;
498  double _eta_cut;
499 };
500 
501 
502 bool Lepjets_Event_Cutter::operator () (const Lepjets_Event_Lep& l) const
503 //
504 // Purpose: Object cut evaluator.
505 //
506 {
507  return ! (l.p().perp() > _pt_cut && abs (l.p().pseudoRapidity()) < _eta_cut);
508 }
509 
510 
511 } // unnamed namespace
512 
513 
514 int Lepjets_Event::cut_leps (double pt_cut, double eta_cut)
515 //
516 // Purpose: Remove all leptons failing the pt and eta cuts.
517 //
518 // Inputs:
519 // pt_cut - Pt cut. Remove objects with pt less than this.
520 // eta_cut - Eta cut. Remove objects with abs(eta) larger than this.
521 //
522 // Returns:
523 // The number of leptons remaining after the cuts.
524 //
525 {
526  _leps.erase (remove_if (_leps.begin(), _leps.end(),
527  Lepjets_Event_Cutter (pt_cut, eta_cut)),
528  _leps.end ());
529  return _leps.size ();
530 }
531 
532 
533 int Lepjets_Event::cut_jets (double pt_cut, double eta_cut)
534 //
535 // Purpose: Remove all jets failing the pt and eta cuts.
536 //
537 // Inputs:
538 // pt_cut - Pt cut. Remove objects with pt less than this.
539 // eta_cut - Eta cut. Remove objects with abs(eta) larger than this.
540 //
541 // Returns:
542 // The number of jets remaining after the cuts.
543 //
544 {
545  _jets.erase (remove_if (_jets.begin(), _jets.end(),
546  Lepjets_Event_Cutter (pt_cut, eta_cut)),
547  _jets.end ());
548  return _jets.size ();
549 }
550 
551 
553 //
554 // Purpose: Remove all but the first N jets.
555 //
556 // Inputs:
557 // n - The number of jets to keep.
558 //
559 {
560  if (n >= _jets.size())
561  return;
562  _jets.erase (_jets.begin() + n, _jets.end());
563 }
564 
565 
566 std::ostream& Lepjets_Event::dump (std::ostream& s, bool full /*=false*/) const
567 //
568 // Purpose: Dump out this object.
569 //
570 // Inputs:
571 // s - The stream to which to write.
572 // full - If true, dump all information for this object.
573 //
574 // Returns:
575 // The stream S.
576 //
577 {
578  s << "Run: " << _runnum << " Event: " << _evnum << "\n";
579  s << "Leptons:\n";
580  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
581  s << " ";
582  _leps[i].dump (s, full);
583  s << "\n";
584  }
585  s << "Jets:\n";
586  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
587  s << " ";
588  _jets[i].dump (s, full);
589  s << "\n";
590  }
591  s << "Missing Et: " << _met << "\n";
592  if (_zvertex != 0)
593  s << "z-vertex: " << _zvertex << "\n";
594  return s;
595 }
596 
597 
599 //
600 // Purpose: Return a string representation of the jet permutation
601 // g - isr/gluon
602 // b - leptonic b
603 // B - hadronic b
604 // w - hadronic W
605 // h - Higgs to b-bbar
606 // ? - Unknown
607 {
608  std::string permutation;
609  for (size_t jet = 0 ; jet != _jets.size() ; ++jet) {
610  permutation = permutation + hitfit::jetTypeString(_jets[jet].type());
611  }
612  return permutation;
613 }
614 
623 std::ostream& operator<< (std::ostream& s, const Lepjets_Event& ev)
624 //
625 // Purpose: Dump out this object.
626 //
627 // Inputs:
628 // s - The stream to which to write.
629 // ev - The object to dump.
630 //
631 // Returns:
632 // The stream S.
633 //
634 {
635  return ev.dump (s);
636 }
637 
638 
639 } // namespace hitfit
640 
641 
642 
643 
644 
const int & runnum() const
Return a constant reference to the run number.
double zvertex() const
Return the value of z-vertex.
type
Definition: HCALResponse.h:21
void trimjets(std::vector< Lepjets_Event_Jet >::size_type n)
Remove all but the first n jets.
Fourvec kt() const
Return the sum of all objects&#39; four-momentum and missing transverse energy.
Represent a lepton in an instance of Lepjets_Event class. This class hold the following information: ...
Lepjets_Event(int runnum, int evnum)
Constructor.
int cut_leps(double pt_cut, double eta_cut)
Remove leptons which fail transverse momentum and pseudorapidity cut.
Calculate and represent resolution for a physical quantity.
Definition: Resolution.h:102
int cut_jets(double pt_cut, double eta_cut)
Remove jets which fail transverse momentum and pseudorapidity cut.
Resolution & kt_res()
Return a reference to the resolution.
bool set_jet_types(const std::vector< int > &)
Set the jet types in the event.
void add_lep(const Lepjets_Event_Lep &lep)
Add a new lepton to the event.
#define X(str)
Definition: MuonsGrabber.cc:48
std::vector< Lepjets_Event_Lep > _leps
bool ev
double pick(double x, double p, CLHEP::HepRandomEngine &engine) const
Generate random value from a Gaussian distribution described by this resolution. Given a value ...
Definition: Resolution.cc:228
uint16_t size_type
Lepjets_Event_Jet & jet(std::vector< Lepjets_Event_Jet >::size_type i)
Return a reference to jet at index position i.
void setMC(bool isMC)
Set the Monte Carlo flag.
double & dnn()
Return a reference to the value of neural network (NN) discriminant (Irrelevant for non-D0 experiment...
int & type()
Return a reference to the type code.
double & dlb()
Return a reference to the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment)...
std::string jetTypeString(int type)
Helper function: Translate jet type code from integer to char. The following notation is used for eac...
Represent a simple event consisting of lepton(s) and jet(s).
std::vector< Lepjets_Event_Jet >::size_type njets() const
Return the number of jets in the event.
std::vector< Lepjets_Event_Lep >::size_type nleps() const
Return the number of leptons in the event.
Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list o...
Definition: Lepjets_Event.h:66
std::vector< int > jet_types() const
Return the jet types in the event.
Definition: GenABIO.cc:168
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:57
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
A class to represent a jet in an instance of Lepjets_Event class. The class is derived from the Lepje...
Fourvec sum(int type) const
Return the sum of all objects&#39; four-momentum which have a particular type.
double _pt_cut
std::ostream & operator<<(std::ostream &s, const Constraint_Intermed &ci)
Output stream operator, print the content of this Constraint_Intermed to an output stream...
Fourvec & met()
Return a reference to the missing transverse energy.
Lepjets_Event_Lep & lep(std::vector< Lepjets_Event_Lep >::size_type i)
Return a reference to lepton at index position i.
std::ostream & dump(std::ostream &s, bool full=false) const
Print the content of this object.
void add_jet(const Lepjets_Event_Jet &jet)
Add a new jet to the event.
bool isMC() const
Return the Monte Carlo flag.
void sort()
Sort objects in the event according to their transverse momentum .
const int & evnum() const
Return a constant reference to the event number.
double _eta_cut
std::vector< Lepjets_Event_Jet > _jets
std::string jet_permutation() const
Return a string representing the jet permutation. The following notation is used for each type of jet...
void smear(CLHEP::HepRandomEngine &engine, bool smear_dir=false)
Smear the objects in the event according to their resolutions.