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 
35 #include <algorithm>
36 #include <functional>
37 #include <cmath>
38 #include <cassert>
39 
40 using std::abs;
41 using std::remove_if;
42 using std::stable_sort;
43 
44 namespace hitfit {
45 
47  //
48  // Purpose: Constructor.
49  //
50  // Inputs:
51  // runnum - The run number.
52  // evnum - The event number.
53  //
54  : _zvertex(0), _isMC(false), _runnum(runnum), _evnum(evnum), _dlb(-1), _dnn(-1) {}
55 
57  //
58  // Purpose: Return the run number.
59  //
60  // Returns:
61  // The run number.
62  //
63  {
64  return _runnum;
65  }
66 
68  //
69  // Purpose: Return the run number.
70  //
71  // Returns:
72  // The run number.
73  //
74  {
75  return _runnum;
76  }
77 
79  //
80  // Purpose: Return the event number.
81  //
82  // Returns:
83  // The event number.
84  //
85  {
86  return _evnum;
87  }
88 
90  //
91  // Purpose: Return the event number.
92  //
93  // Returns:
94  // The event number.
95  //
96  {
97  return _evnum;
98  }
99 
101  //
102  // Purpose: Return the length of the lepton list.
103  //
104  // Returns:
105  // The length of the lepton list.
106  //
107  {
108  return _leps.size();
109  }
110 
112  //
113  // Purpose: Return the length of the jet list.
114  //
115  // Returns:
116  // The length of the jet list.
117  //
118  {
119  return _jets.size();
120  }
121 
123  //
124  // Purpose: Return the Ith lepton.
125  //
126  // Inputs:
127  // i - The lepton index (counting from 0).
128  //
129  // Returns:
130  // The Ith lepton.
131  //
132  {
133  assert(i < _leps.size());
134  return _leps[i];
135  }
136 
138  //
139  // Purpose: Return the Ith jet.
140  //
141  // Inputs:
142  // i - The jet index (counting from 0).
143  //
144  // Returns:
145  // The Ith jet.
146  //
147  {
148  assert(i < _jets.size());
149  return _jets[i];
150  }
151 
153  //
154  // Purpose: Return the Ith lepton.
155  //
156  // Inputs:
157  // i - The lepton index (counting from 0).
158  //
159  // Returns:
160  // The Ith lepton.
161  //
162  {
163  assert(i < _leps.size());
164  return _leps[i];
165  }
166 
168  //
169  // Purpose: Return the Ith jet.
170  //
171  // Inputs:
172  // i - The jet index (counting from 0).
173  //
174  // Returns:
175  // The Ith jet.
176  //
177  {
178  assert(i < _jets.size());
179  return _jets[i];
180  }
181 
183  //
184  // Purpose: Return the missing Et.
185  //
186  // Returns:
187  // The missing Et.
188  //
189  {
190  return _met;
191  }
192 
194  //
195  // Purpose: Return the missing Et.
196  //
197  // Returns:
198  // The missing Et.
199  //
200  {
201  return _met;
202  }
203 
205  //
206  // Purpose: Return the kt resolution.
207  //
208  // Returns:
209  // The kt resolution.
210  //
211  {
212  return _kt_res;
213  }
214 
216  //
217  // Purpose: Return the kt resolution.
218  //
219  // Returns:
220  // The kt resolution.
221  //
222  {
223  return _kt_res;
224  }
225 
227  //
228  // Purpose: Return the z-vertex.
229  //
230  // Returns:
231  // The z-vertex.
232  //
233  {
234  return _zvertex;
235  }
236 
238  //
239  // Purpose: Return the z-vertex.
240  //
241  // Returns:
242  // The z-vertex.
243  //
244  {
245  return _zvertex;
246  }
247 
249  //
250  // Purpose: Return the isMC flag.
251  //
252  // Returns:
253  // The isMC flag.
254  //
255  {
256  return _isMC;
257  }
258 
260  //
261  // Purpose: set isMC flag.
262  //
263  // Returns:
264  // nothing
265  //
266  {
267  _isMC = isMC;
268  }
269 
271  //
272  // Purpose: Return the LB discriminant.
273  //
274  // Returns:
275  // The LB discriminant.
276  //
277  {
278  return _dlb;
279  }
280 
282  //
283  // Purpose: Return the LB discriminant.
284  //
285  // Returns:
286  // The LB discriminant.
287  //
288  {
289  return _dlb;
290  }
291 
293  //
294  // Purpose: Return the NN discriminant.
295  //
296  // Returns:
297  // The NN discriminant.
298  //
299  {
300  return _dnn;
301  }
302 
304  //
305  // Purpose: Return the NN discriminant.
306  //
307  // Returns:
308  // The NN discriminant.
309  //
310  {
311  return _dnn;
312  }
313 
315  //
316  // Purpose: Sum all objects with type code TYPE.
317  //
318  // Inputs:
319  // type - The type code to match.
320  //
321  // Returns:
322  // The sum of all objects with type code TYPE.
323  //
324  {
325  Fourvec out;
326  for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++)
327  if (_leps[i].type() == type)
328  out += _leps[i].p();
329  for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++)
330  if (_jets[i].type() == type)
331  out += _jets[i].p();
332  return out;
333  }
334 
336  //
337  // Purpose: Calculate kt --- sum of all objects plus missing Et.
338  //
339  // Returns:
340  // The event kt.
341  {
342  Fourvec v = _met;
343  for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++)
344  v += _leps[i].p();
345  for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++)
346  v += _jets[i].p();
347  return v;
348  }
349 
351  //
352  // Purpose: Add a lepton to the event.
353  //
354  // Inputs:
355  // lep - The lepton to add.
356  //
357  {
358  _leps.push_back(lep);
359  }
360 
362  //
363  // Purpose: Add a jet to the event.
364  //
365  // Inputs:
366  // jet - The jet to add.
367  //
368  {
369  _jets.push_back(jet);
370  }
371 
372  void Lepjets_Event::smear(CLHEP::HepRandomEngine& engine, bool smear_dir /*= false*/)
373  //
374  // Purpose: Smear the objects in the event according to their resolutions.
375  //
376  // Inputs:
377  // engine - The underlying RNG.
378  // smear_dir - If false, smear the momentum only.
379  //
380  {
381  Fourvec before, after;
382  for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++) {
383  before += _leps[i].p();
384  _leps[i].smear(engine, smear_dir);
385  after += _leps[i].p();
386  }
387  for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++) {
388  before += _jets[i].p();
389  _jets[i].smear(engine, smear_dir);
390  after += _jets[i].p();
391  }
392 
393  Fourvec kt = _met + before;
396  _met = kt - after;
397  }
398 
400  //
401  // Purpose: Sort the objects in the event in order of descending pt.
402  //
403  {
404  std::stable_sort(_leps.begin(), _leps.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
405  std::stable_sort(_jets.begin(), _jets.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
406  }
407 
408  std::vector<int> Lepjets_Event::jet_types() const
409  //
410  // Purpose: Return the jet types of the event
411  //
412  {
413  std::vector<int> ret;
414  for (std::vector<Lepjets_Event_Jet>::size_type ijet = 0; ijet != _jets.size(); ijet++) {
415  ret.push_back(jet(ijet).type());
416  }
417  return ret;
418  }
419 
420  bool Lepjets_Event::set_jet_types(const std::vector<int>& _jet_types)
421  //
422  // Purpose: Set the jet types of the event
423  // Return false if it fails, trus if it succeeds
424  //
425  {
426  if (_jets.size() != _jet_types.size()) {
427  return false;
428  }
429  bool saw_hadw1 = false;
431  int t = _jet_types[i];
432  if (t == hadw1_label) {
433  if (saw_hadw1)
434  t = hadw2_label;
435  saw_hadw1 = true;
436  }
437  jet(i).type() = t;
438  }
439  return true;
440  }
441 
442  namespace {
443 
444  struct Lepjets_Event_Cutter
445  //
446  // Purpose: Helper for cutting on objects.
447  //
448  {
449  Lepjets_Event_Cutter(double pt_cut, double eta_cut) : _pt_cut(pt_cut), _eta_cut(eta_cut) {}
450  bool operator()(const Lepjets_Event_Lep& l) const;
451  double _pt_cut;
452  double _eta_cut;
453  };
454 
455  bool Lepjets_Event_Cutter::operator()(const Lepjets_Event_Lep& l) const
456  //
457  // Purpose: Object cut evaluator.
458  //
459  {
460  return !(l.p().perp() > _pt_cut && abs(l.p().pseudoRapidity()) < _eta_cut);
461  }
462 
463  } // unnamed namespace
464 
465  int Lepjets_Event::cut_leps(double pt_cut, double eta_cut)
466  //
467  // Purpose: Remove all leptons failing the pt and eta cuts.
468  //
469  // Inputs:
470  // pt_cut - Pt cut. Remove objects with pt less than this.
471  // eta_cut - Eta cut. Remove objects with abs(eta) larger than this.
472  //
473  // Returns:
474  // The number of leptons remaining after the cuts.
475  //
476  {
477  _leps.erase(remove_if(_leps.begin(), _leps.end(), Lepjets_Event_Cutter(pt_cut, eta_cut)), _leps.end());
478  return _leps.size();
479  }
480 
481  int Lepjets_Event::cut_jets(double pt_cut, double eta_cut)
482  //
483  // Purpose: Remove all jets failing the pt and eta cuts.
484  //
485  // Inputs:
486  // pt_cut - Pt cut. Remove objects with pt less than this.
487  // eta_cut - Eta cut. Remove objects with abs(eta) larger than this.
488  //
489  // Returns:
490  // The number of jets remaining after the cuts.
491  //
492  {
493  _jets.erase(remove_if(_jets.begin(), _jets.end(), Lepjets_Event_Cutter(pt_cut, eta_cut)), _jets.end());
494  return _jets.size();
495  }
496 
498  //
499  // Purpose: Remove all but the first N jets.
500  //
501  // Inputs:
502  // n - The number of jets to keep.
503  //
504  {
505  if (n >= _jets.size())
506  return;
507  _jets.erase(_jets.begin() + n, _jets.end());
508  }
509 
510  std::ostream& Lepjets_Event::dump(std::ostream& s, bool full /*=false*/) const
511  //
512  // Purpose: Dump out this object.
513  //
514  // Inputs:
515  // s - The stream to which to write.
516  // full - If true, dump all information for this object.
517  //
518  // Returns:
519  // The stream S.
520  //
521  {
522  s << "Run: " << _runnum << " Event: " << _evnum << "\n";
523  s << "Leptons:\n";
524  for (std::vector<Lepjets_Event_Lep>::size_type i = 0; i < _leps.size(); i++) {
525  s << " ";
526  _leps[i].dump(s, full);
527  s << "\n";
528  }
529  s << "Jets:\n";
530  for (std::vector<Lepjets_Event_Jet>::size_type i = 0; i < _jets.size(); i++) {
531  s << " ";
532  _jets[i].dump(s, full);
533  s << "\n";
534  }
535  s << "Missing Et: " << _met << "\n";
536  if (_zvertex != 0)
537  s << "z-vertex: " << _zvertex << "\n";
538  return s;
539  }
540 
542  //
543  // Purpose: Return a string representation of the jet permutation
544  // g - isr/gluon
545  // b - leptonic b
546  // B - hadronic b
547  // w - hadronic W
548  // h - Higgs to b-bbar
549  // ? - Unknown
550  {
551  std::string permutation;
552  for (size_t jet = 0; jet != _jets.size(); ++jet) {
553  permutation = permutation + hitfit::jetTypeString(_jets[jet].type());
554  }
555  return permutation;
556  }
557 
566  std::ostream& operator<<(std::ostream& s, const Lepjets_Event& ev)
567  //
568  // Purpose: Dump out this object.
569  //
570  // Inputs:
571  // s - The stream to which to write.
572  // ev - The object to dump.
573  //
574  // Returns:
575  // The stream S.
576  //
577  {
578  return ev.dump(s);
579  }
580 
581 } // namespace hitfit
void trimjets(std::vector< Lepjets_Event_Jet >::size_type n)
Remove all but the first n jets.
double zvertex() const
Return the value of z-vertex.
int runnum
Represent a lepton in an instance of Lepjets_Event class. This class hold the following information: ...
Lepjets_Event(int runnum, int evnum)
Constructor.
const int & evnum() const
Return a constant reference to the event number.
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:182
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:98
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.
ret
prodAgent to be discontinued
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:38
std::vector< Lepjets_Event_Lep >::size_type nleps() const
Return the number of leptons in the event.
std::vector< Lepjets_Event_Lep > _leps
assert(be >=bs)
std::ostream & dump(std::ostream &s, bool full=false) const
Print the content of this object.
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...
bool isMC() const
Return the Monte Carlo flag.
int & type()
Return a reference to the type code.
Fourvec sum(int type) const
Return the sum of all objects&#39; four-momentum which have a particular type.
double & dlb()
Return a reference to the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment)...
dictionary isMC
Definition: PV_cfg.py:29
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).
Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list o...
Definition: Lepjets_Event.h:62
Definition: GenABIO.cc:168
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:55
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...
const int & runnum() const
Return a constant reference to the run number.
double _pt_cut
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::string jet_permutation() const
Return a string representing the jet permutation. The following notation is used for each type of jet...
std::vector< int > jet_types() const
Return the jet types in the event.
void add_jet(const Lepjets_Event_Jet &jet)
Add a new jet to the event.
void sort()
Sort objects in the event according to their transverse momentum .
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 kt() const
Return the sum of all objects&#39; four-momentum and missing transverse energy.
double _eta_cut
std::vector< Lepjets_Event_Jet > _jets
void smear(CLHEP::HepRandomEngine &engine, bool smear_dir=false)
Smear the objects in the event according to their resolutions.
std::vector< Lepjets_Event_Jet >::size_type njets() const
Return the number of jets in the event.