CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Lepjets_Event.h
Go to the documentation of this file.
1 //
2 //
3 // File: hitfit/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 // An instance of this class holds a list of `leptons' (as represented
8 // by Lepjets_Event_Lep) and `jets' (as represented by Lepjets_Event_Jet).
9 // We also record:
10 //
11 // - Missing Et
12 // - z-vertex
13 // - run and event number
14 //
15 // CMSSW File : interface/Lepjets_Event.h
16 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
17 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
18 //
19 
41 #ifndef HITFIT_LEPJETS_EVENT_H
42 #define HITFIT_LEPJETS_EVENT_H
43 
46 #include <vector>
47 #include <iosfwd>
48 
49 namespace hitfit {
50 
63  //
64  // Purpose: Represent a simple `event' consisting of leptons and jets.
65  //
66  {
67  public:
68  // Constructor.
76  Lepjets_Event(int runnum, int evnum);
77 
78  // Get the run and event number.
79 
83  const int& runnum() const;
84 
88  int& runnum();
89 
93  const int& evnum() const;
94 
98  int& evnum();
99 
100  // Get the length of the lepton and jet lists.
105 
110 
111  // Access leptons and jets.
112 
119 
126 
133 
140 
141  // Access missing Et.
145  Fourvec& met();
146 
150  const Fourvec& met() const;
151 
152  // Access kt resolution.
153 
157  Resolution& kt_res();
158 
162  const Resolution& kt_res() const;
163 
164  // Access the z-vertex.
165 
169  double zvertex() const;
170 
174  double& zvertex();
175 
176  // Access the isMC flag.
180  bool isMC() const;
181 
185  void setMC(bool isMC);
186 
187  // Access the discriminants.
192  double& dlb();
193 
198  double dlb() const;
199 
204  double& dnn();
205 
210  double dnn() const;
211 
212  // Sum all objects (leptons or jets) with type TYPE.
219  Fourvec sum(int type) const;
220 
221  // Calculate kt --- sum of all objects plus missing Et.
226  Fourvec kt() const;
227 
228  // Add new objects to the event.
234  void add_lep(const Lepjets_Event_Lep& lep);
235 
241  void add_jet(const Lepjets_Event_Jet& jet);
242 
243  // Smear the objects in the event according to their resolutions.
252  void smear(CLHEP::HepRandomEngine& engine, bool smear_dir = false);
253 
254  // Sort according to pt.
259  void sort();
260 
261  // Get jet types
265  std::vector<int> jet_types() const;
266 
267  // Set jet types
271  bool set_jet_types(const std::vector<int>&);
272 
273  // Remove objects failing pt and eta cuts.
284  int cut_leps(double pt_cut, double eta_cut);
285 
296  int cut_jets(double pt_cut, double eta_cut);
297 
298  // Remove all but the first N jets.
305 
306  // Dump this object.
317  std::ostream& dump(std::ostream& s, bool full = false) const;
318 
330 
331  private:
332  // The lepton and jet lists.
333 
337  std::vector<Lepjets_Event_Lep> _leps;
338 
342  std::vector<Lepjets_Event_Jet> _jets;
343 
344  // Other event state.
349 
354 
358  double _zvertex;
359 
363  bool _isMC;
364 
368  int _runnum;
369 
373  int _evnum;
374 
378  double _dlb;
379 
383  double _dnn;
384  };
385 
386  // Print the object.
387  std::ostream& operator<<(std::ostream& s, const Lepjets_Event& ev);
388 
389 } // namespace hitfit
390 
391 #endif // not HITFIT_LEPJETS_EVENT_H
const int & runnum() const
Return a constant reference to the run number.
double zvertex() const
Return the value of z-vertex.
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: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.
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.
std::vector< Lepjets_Event_Lep > _leps
bool ev
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...
double & dlb()
Return a reference to the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment)...
Represent a jet in an instance of Lepjets_Event class.
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:62
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:55
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.
Represent a lepton in an instance of Lepjets_Event class.
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.
std::ostream & operator<<(std::ostream &s, const Constraint_Intermed &ci)
Output stream operator, print the content of this Constraint_Intermed to an output stream...
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.