CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Lepjets_Event.h
Go to the documentation of this file.
1 //
2 // $Id: Lepjets_Event.h,v 1.1 2011/05/26 09:46:53 mseidel Exp $
3 //
4 // File: hitfit/Lepjets_Event.h
5 // Purpose: Represent a simple `event' consisting of leptons and jets.
6 // Created: Jul, 2000, sss, based on run 1 mass analysis code.
7 //
8 // An instance of this class holds a list of `leptons' (as represented
9 // by Lepjets_Event_Lep) and `jets' (as represented by Lepjets_Event_Jet).
10 // We also record:
11 //
12 // - Missing Et
13 // - z-vertex
14 // - run and event number
15 //
16 // CMSSW File : interface/Lepjets_Event.h
17 // Original Author : Scott Stuart Snyder <snyder@bnl.gov> for D0
18 // Imported to CMSSW by Haryo Sumowidagdo <Suharyo.Sumowidagdo@cern.ch>
19 //
20 
21 
43 #ifndef HITFIT_LEPJETS_EVENT_H
44 #define HITFIT_LEPJETS_EVENT_H
45 
46 
49 #include <vector>
50 #include <iosfwd>
51 
52 
53 namespace hitfit {
54 
55 
68 //
69 // Purpose: Represent a simple `event' consisting of leptons and jets.
70 //
71 {
72 public:
73 
74  // Constructor.
82  Lepjets_Event (int runnum, int evnum);
83 
84  // Get the run and event number.
85 
89  const int& runnum () const;
90 
94  int& runnum();
95 
99  const int& evnum () const;
100 
104  int& evnum();
105 
106  // Get the length of the lepton and jet lists.
111 
116 
117  // Access leptons and jets.
118 
125 
132 
139 
146 
147  // Access missing Et.
151  Fourvec& met ();
152 
156  const Fourvec& met () const;
157 
158  // Access kt resolution.
159 
163  Resolution& kt_res ();
164 
168  const Resolution& kt_res () const;
169 
170  // Access the z-vertex.
171 
175  double zvertex () const;
176 
180  double& zvertex ();
181 
182  // Access the isMC flag.
186  bool isMC () const;
187 
191  void setMC (bool isMC);
192 
193  // Access the discriminants.
198  double& dlb ();
199 
204  double dlb () const;
205 
210  double& dnn ();
211 
216  double dnn () const;
217 
218  // Sum all objects (leptons or jets) with type TYPE.
225  Fourvec sum (int type) const;
226 
227  // Calculate kt --- sum of all objects plus missing Et.
232  Fourvec kt () const;
233 
234  // Add new objects to the event.
240  void add_lep (const Lepjets_Event_Lep& lep);
241 
247  void add_jet (const Lepjets_Event_Jet& jet);
248 
249  // Smear the objects in the event according to their resolutions.
258  void smear (CLHEP::HepRandomEngine& engine, bool smear_dir = false);
259 
260  // Sort according to pt.
265  void sort ();
266 
267  // Get jet types
271  std::vector<int> jet_types() const;
272 
273  // Set jet types
277  bool set_jet_types(const std::vector<int>&);
278 
279  // Remove objects failing pt and eta cuts.
290  int cut_leps (double pt_cut, double eta_cut);
291 
302  int cut_jets (double pt_cut, double eta_cut);
303 
304  // Remove all but the first N jets.
311 
312  // Dump this object.
323  std::ostream& dump (std::ostream& s, bool full = false) const;
324 
336 
337 private:
338  // The lepton and jet lists.
339 
343  std::vector<Lepjets_Event_Lep> _leps;
344 
348  std::vector<Lepjets_Event_Jet> _jets;
349 
350  // Other event state.
355 
360 
364  double _zvertex;
365 
369  bool _isMC;
370 
374  int _runnum;
375 
379  int _evnum;
380 
384  double _dlb;
385 
389  double _dnn;
390 };
391 
392 
393 // Print the object.
394 std::ostream& operator<< (std::ostream& s, const Lepjets_Event& ev);
395 
396 
397 
398 } // namespace hitfit
399 
400 
401 #endif // not HITFIT_LEPJETS_EVENT_H
402 
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
int i
Definition: DBlmapReader.cc:9
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:103
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
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:67
std::vector< int > jet_types() const
Return the jet types in the event.
Definition: GenABIO.cc:193
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:58
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.