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 //
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 
20 
42 #ifndef HITFIT_LEPJETS_EVENT_H
43 #define HITFIT_LEPJETS_EVENT_H
44 
45 
48 #include <vector>
49 #include <iosfwd>
50 
51 
52 namespace hitfit {
53 
54 
67 //
68 // Purpose: Represent a simple `event' consisting of leptons and jets.
69 //
70 {
71 public:
72 
73  // Constructor.
81  Lepjets_Event (int runnum, int evnum);
82 
83  // Get the run and event number.
84 
88  const int& runnum () const;
89 
93  int& runnum();
94 
98  const int& evnum () const;
99 
103  int& evnum();
104 
105  // Get the length of the lepton and jet lists.
110 
115 
116  // Access leptons and jets.
117 
124 
131 
138 
145 
146  // Access missing Et.
150  Fourvec& met ();
151 
155  const Fourvec& met () const;
156 
157  // Access kt resolution.
158 
162  Resolution& kt_res ();
163 
167  const Resolution& kt_res () const;
168 
169  // Access the z-vertex.
170 
174  double zvertex () const;
175 
179  double& zvertex ();
180 
181  // Access the isMC flag.
185  bool isMC () const;
186 
190  void setMC (bool isMC);
191 
192  // Access the discriminants.
197  double& dlb ();
198 
203  double dlb () const;
204 
209  double& dnn ();
210 
215  double dnn () const;
216 
217  // Sum all objects (leptons or jets) with type TYPE.
224  Fourvec sum (int type) const;
225 
226  // Calculate kt --- sum of all objects plus missing Et.
231  Fourvec kt () const;
232 
233  // Add new objects to the event.
239  void add_lep (const Lepjets_Event_Lep& lep);
240 
246  void add_jet (const Lepjets_Event_Jet& jet);
247 
248  // Smear the objects in the event according to their resolutions.
257  void smear (CLHEP::HepRandomEngine& engine, bool smear_dir = false);
258 
259  // Sort according to pt.
264  void sort ();
265 
266  // Get jet types
270  std::vector<int> jet_types() const;
271 
272  // Set jet types
276  bool set_jet_types(const std::vector<int>&);
277 
278  // Remove objects failing pt and eta cuts.
289  int cut_leps (double pt_cut, double eta_cut);
290 
301  int cut_jets (double pt_cut, double eta_cut);
302 
303  // Remove all but the first N jets.
310 
311  // Dump this object.
322  std::ostream& dump (std::ostream& s, bool full = false) const;
323 
335 
336 private:
337  // The lepton and jet lists.
338 
342  std::vector<Lepjets_Event_Lep> _leps;
343 
347  std::vector<Lepjets_Event_Jet> _jets;
348 
349  // Other event state.
354 
359 
363  double _zvertex;
364 
368  bool _isMC;
369 
373  int _runnum;
374 
378  int _evnum;
379 
383  double _dlb;
384 
388  double _dnn;
389 };
390 
391 
392 // Print the object.
393 std::ostream& operator<< (std::ostream& s, const Lepjets_Event& ev);
394 
395 
396 
397 } // namespace hitfit
398 
399 
400 #endif // not HITFIT_LEPJETS_EVENT_H
401 
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: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.
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:66
std::vector< int > jet_types() const
Return the jet types in the event.
Definition: GenABIO.cc:180
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:57
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.