test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
hitfit::Lepjets_Event Class Reference

Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list of leptons (as represented by the Lepjets_Event_Lep class) and jets (as represented by Lepjets_Event_Jet class). Also recorded are: More...

#include <Lepjets_Event.h>

Public Member Functions

void add_jet (const Lepjets_Event_Jet &jet)
 Add a new jet to the event. More...
 
void add_lep (const Lepjets_Event_Lep &lep)
 Add a new lepton to the event. More...
 
int cut_jets (double pt_cut, double eta_cut)
 Remove jets which fail transverse momentum $ p_{T} $ and pseudorapidity $ \eta $ cut. More...
 
int cut_leps (double pt_cut, double eta_cut)
 Remove leptons which fail transverse momentum $ p_{T} $ and pseudorapidity $ \eta $ cut. More...
 
double & dlb ()
 Return a reference to the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment). More...
 
double dlb () const
 Return the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment). More...
 
double & dnn ()
 Return a reference to the value of neural network (NN) discriminant (Irrelevant for non-D0 experiment). More...
 
double dnn () const
 Return a the value of neural network (NN) discriminant (Irrelevant for non-D0 experiment). More...
 
std::ostream & dump (std::ostream &s, bool full=false) const
 Print the content of this object. More...
 
const int & evnum () const
 Return a constant reference to the event number. More...
 
int & evnum ()
 Return a reference to the event number. More...
 
bool isMC () const
 Return the Monte Carlo flag. More...
 
Lepjets_Event_Jetjet (std::vector< Lepjets_Event_Jet >::size_type i)
 Return a reference to jet at index position i. More...
 
const Lepjets_Event_Jetjet (std::vector< Lepjets_Event_Jet >::size_type i) const
 Return a constant reference to jet at index position i. More...
 
std::string jet_permutation () const
 Return a string representing the jet permutation. The following notation is used for each type of jet: More...
 
std::vector< int > jet_types () const
 Return the jet types in the event. More...
 
Fourvec kt () const
 Return the sum of all objects' four-momentum and missing transverse energy. More...
 
Resolutionkt_res ()
 Return a reference to the $ k_{T} $ resolution. More...
 
const Resolutionkt_res () const
 Return a const reference to the $ k_{T} $ resolution. More...
 
Lepjets_Event_Leplep (std::vector< Lepjets_Event_Lep >::size_type i)
 Return a reference to lepton at index position i. More...
 
const Lepjets_Event_Leplep (std::vector< Lepjets_Event_Lep >::size_type i) const
 Return a constant reference to lepton at index position i. More...
 
 Lepjets_Event (int runnum, int evnum)
 Constructor. More...
 
Fourvecmet ()
 Return a reference to the missing transverse energy. More...
 
const Fourvecmet () const
 Return a constant reference to the missing transverse energy. More...
 
std::vector< Lepjets_Event_Jet >
::size_type 
njets () const
 Return the number of jets in the event. More...
 
std::vector< Lepjets_Event_Lep >
::size_type 
nleps () const
 Return the number of leptons in the event. More...
 
const int & runnum () const
 Return a constant reference to the run number. More...
 
int & runnum ()
 Return a reference to the run number. More...
 
bool set_jet_types (const std::vector< int > &)
 Set the jet types in the event. More...
 
void setMC (bool isMC)
 Set the Monte Carlo flag. More...
 
void smear (CLHEP::HepRandomEngine &engine, bool smear_dir=false)
 Smear the objects in the event according to their resolutions. More...
 
void sort ()
 Sort objects in the event according to their transverse momentum $ p_{T} $ . More...
 
Fourvec sum (int type) const
 Return the sum of all objects' four-momentum which have a particular type. More...
 
void trimjets (std::vector< Lepjets_Event_Jet >::size_type n)
 Remove all but the first n jets. More...
 
double zvertex () const
 Return the value of z-vertex. More...
 
double & zvertex ()
 Return a reference to the value of z-vertex. More...
 

Private Attributes

double _dlb
 
double _dnn
 
int _evnum
 
bool _isMC
 
std::vector< Lepjets_Event_Jet_jets
 
Resolution _kt_res
 
std::vector< Lepjets_Event_Lep_leps
 
Fourvec _met
 
int _runnum
 
double _zvertex
 

Detailed Description

Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list of leptons (as represented by the Lepjets_Event_Lep class) and jets (as represented by Lepjets_Event_Jet class). Also recorded are:

Definition at line 66 of file Lepjets_Event.h.

Constructor & Destructor Documentation

hitfit::Lepjets_Event::Lepjets_Event ( int  runnum,
int  evnum 
)

Constructor.

Parameters
runnumThe run number.
evnumThe event number.

Definition at line 52 of file Lepjets_Event.cc.

60  : _zvertex (0),
61  _isMC(false),
62  _runnum (runnum),
63  _evnum (evnum),
64  _dlb (-1),
65  _dnn (-1)
66 {
67 }
const int & runnum() const
Return a constant reference to the run number.
const int & evnum() const
Return a constant reference to the event number.

Member Function Documentation

void hitfit::Lepjets_Event::add_jet ( const Lepjets_Event_Jet jet)

Add a new jet to the event.

Parameters
jetThe jet to be added.

Definition at line 399 of file Lepjets_Event.cc.

References _jets.

Referenced by hitfit::RunHitFit< AElectron, AMuon, AJet, AMet >::FitAllPermutation(), hitfit::gentop(), and hitfit::gentth().

406 {
407  _jets.push_back (jet);
408 }
Lepjets_Event_Jet & jet(std::vector< Lepjets_Event_Jet >::size_type i)
Return a reference to jet at index position i.
std::vector< Lepjets_Event_Jet > _jets
void hitfit::Lepjets_Event::add_lep ( const Lepjets_Event_Lep lep)

Add a new lepton to the event.

Parameters
lepThe lepton to be added.

Definition at line 387 of file Lepjets_Event.cc.

References _leps.

Referenced by hitfit::RunHitFit< AElectron, AMuon, AJet, AMet >::AddLepton(), hitfit::gentop(), and hitfit::gentth().

394 {
395  _leps.push_back (lep);
396 }
std::vector< Lepjets_Event_Lep > _leps
Lepjets_Event_Lep & lep(std::vector< Lepjets_Event_Lep >::size_type i)
Return a reference to lepton at index position i.
int hitfit::Lepjets_Event::cut_jets ( double  pt_cut,
double  eta_cut 
)

Remove jets which fail transverse momentum $ p_{T} $ and pseudorapidity $ \eta $ cut.

Parameters
pt_cutRemove jets which have transverse momentum $ p_{T} $ less than this value, in GeV.
eta_cutRemove jetss which have absolute pseudorapidity $ |\eta| $ more than this value.

Definition at line 535 of file Lepjets_Event.cc.

References _jets.

546 {
547  _jets.erase (remove_if (_jets.begin(), _jets.end(),
548  Lepjets_Event_Cutter (pt_cut, eta_cut)),
549  _jets.end ());
550  return _jets.size ();
551 }
std::vector< Lepjets_Event_Jet > _jets
int hitfit::Lepjets_Event::cut_leps ( double  pt_cut,
double  eta_cut 
)

Remove leptons which fail transverse momentum $ p_{T} $ and pseudorapidity $ \eta $ cut.

Parameters
pt_cutRemove leptons which have transverse momentum $ p_{T} $ less than this value, in GeV.
eta_cutRemove leptons which have absolute pseudorapidity $ |\eta| $ more than this value.

Definition at line 516 of file Lepjets_Event.cc.

References _leps.

527 {
528  _leps.erase (remove_if (_leps.begin(), _leps.end(),
529  Lepjets_Event_Cutter (pt_cut, eta_cut)),
530  _leps.end ());
531  return _leps.size ();
532 }
std::vector< Lepjets_Event_Lep > _leps
double & hitfit::Lepjets_Event::dlb ( )

Return a reference to the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment).

Definition at line 313 of file Lepjets_Event.cc.

References _dlb.

320 {
321  return _dlb;
322 }
double hitfit::Lepjets_Event::dlb ( ) const

Return the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment).

Definition at line 301 of file Lepjets_Event.cc.

References _dlb.

308 {
309  return _dlb;
310 }
double & hitfit::Lepjets_Event::dnn ( )

Return a reference to the value of neural network (NN) discriminant (Irrelevant for non-D0 experiment).

Definition at line 337 of file Lepjets_Event.cc.

References _dnn.

344 {
345  return _dnn;
346 }
double hitfit::Lepjets_Event::dnn ( ) const

Return a the value of neural network (NN) discriminant (Irrelevant for non-D0 experiment).

Definition at line 325 of file Lepjets_Event.cc.

References _dnn.

332 {
333  return _dnn;
334 }
std::ostream & hitfit::Lepjets_Event::dump ( std::ostream &  s,
bool  full = false 
) const

Print the content of this object.

Parameters
sThe output stream to which to write
fullIf TRUE, print all information about this instance of Lepjets_Event.
If FALSE, print partial information about this instance of Lepjets_Event.

Definition at line 568 of file Lepjets_Event.cc.

References full, i, and alignCSCRings::s.

Referenced by hitfit::operator<<().

579 {
580  s << "Run: " << _runnum << " Event: " << _evnum << "\n";
581  s << "Leptons:\n";
582  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
583  s << " ";
584  _leps[i].dump (s, full);
585  s << "\n";
586  }
587  s << "Jets:\n";
588  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
589  s << " ";
590  _jets[i].dump (s, full);
591  s << "\n";
592  }
593  s << "Missing Et: " << _met << "\n";
594  if (_zvertex != 0)
595  s << "z-vertex: " << _zvertex << "\n";
596  return s;
597 }
int i
Definition: DBlmapReader.cc:9
std::vector< Lepjets_Event_Lep > _leps
uint16_t size_type
Definition: GenABIO.cc:180
std::vector< Lepjets_Event_Jet > _jets
const int & hitfit::Lepjets_Event::evnum ( ) const

Return a constant reference to the event number.

Definition at line 94 of file Lepjets_Event.cc.

References _evnum.

101 {
102  return _evnum;
103 }
int & hitfit::Lepjets_Event::evnum ( )

Return a reference to the event number.

Definition at line 106 of file Lepjets_Event.cc.

References _evnum.

113 {
114  return _evnum;
115 }
bool hitfit::Lepjets_Event::isMC ( ) const

Return the Monte Carlo flag.

Definition at line 278 of file Lepjets_Event.cc.

References _isMC.

Referenced by setMC().

285 {
286  return _isMC;
287 }
Lepjets_Event_Jet & hitfit::Lepjets_Event::jet ( std::vector< Lepjets_Event_Jet >::size_type  i)

Return a reference to jet at index position i.

Parameters
iThe jet index position.

Definition at line 158 of file Lepjets_Event.cc.

References _jets, assert(), and i.

Referenced by jet_permutation(), jet_types(), TtSemiLepHitFitProducer< LeptonCollection >::produce(), and set_jet_types().

168 {
169  assert (i < _jets.size());
170  return _jets[i];
171 }
int i
Definition: DBlmapReader.cc:9
assert(m_qm.get())
std::vector< Lepjets_Event_Jet > _jets
const Lepjets_Event_Jet & hitfit::Lepjets_Event::jet ( std::vector< Lepjets_Event_Jet >::size_type  i) const

Return a constant reference to jet at index position i.

Parameters
iThe jet index position.

Definition at line 190 of file Lepjets_Event.cc.

References assert(), and i.

200 {
201  assert (i < _jets.size());
202  return _jets[i];
203 }
int i
Definition: DBlmapReader.cc:9
assert(m_qm.get())
std::vector< Lepjets_Event_Jet > _jets
std::string hitfit::Lepjets_Event::jet_permutation ( ) const

Return a string representing the jet permutation. The following notation is used for each type of jet:

  • g ISR/gluon.
  • b leptonic $ b- $ quark.
  • B hadronic $ b- $ quark.
  • w hadronic jet from $ W- $ boson.
  • H $ b- $ jet from Higgs boson.
  • ? Unknown.

Definition at line 600 of file Lepjets_Event.cc.

References _jets, jet(), hitfit::jetTypeString(), and AlCaHLTBitMon_QueryRunRegistry::string.

609 {
610  std::string permutation;
611  for (size_t jet = 0 ; jet != _jets.size() ; ++jet) {
612  permutation = permutation + hitfit::jetTypeString(_jets[jet].type());
613  }
614  return permutation;
615 }
type
Definition: HCALResponse.h:21
Lepjets_Event_Jet & jet(std::vector< Lepjets_Event_Jet >::size_type i)
Return a reference to jet at index position i.
std::string jetTypeString(int type)
Helper function: Translate jet type code from integer to char. The following notation is used for eac...
std::vector< Lepjets_Event_Jet > _jets
std::vector< int > hitfit::Lepjets_Event::jet_types ( ) const

Return the jet types in the event.

Definition at line 449 of file Lepjets_Event.cc.

References _jets, jet(), and run_regression::ret.

Referenced by hitfit::Fit_Result::jet_types().

453 {
454  std::vector<int> ret;
456  ijet != _jets.size() ;
457  ijet++) {
458  ret.push_back(jet(ijet).type());
459  }
460  return ret;
461 }
type
Definition: HCALResponse.h:21
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.
std::vector< Lepjets_Event_Jet > _jets
Fourvec hitfit::Lepjets_Event::kt ( ) const

Return the sum of all objects' four-momentum and missing transverse energy.

Definition at line 371 of file Lepjets_Event.cc.

References _jets, _leps, _met, i, AlCaHLTBitMon_ParallelJobs::p, and findQualityFiles::v.

Referenced by smear().

377 {
378  Fourvec v = _met;
380  v += _leps[i].p();
382  v += _jets[i].p();
383  return v;
384 }
int i
Definition: DBlmapReader.cc:9
std::vector< Lepjets_Event_Lep > _leps
uint16_t size_type
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:57
std::vector< Lepjets_Event_Jet > _jets
Resolution & hitfit::Lepjets_Event::kt_res ( )

Return a reference to the $ k_{T} $ resolution.

Definition at line 230 of file Lepjets_Event.cc.

References _kt_res.

Referenced by hitfit::gentop(), hitfit::gentth(), hitfit::RunHitFit< AElectron, AMuon, AJet, AMet >::SetKtResolution(), and hitfit::RunHitFit< AElectron, AMuon, AJet, AMet >::SetMet().

237 {
238  return _kt_res;
239 }
const Resolution & hitfit::Lepjets_Event::kt_res ( ) const

Return a const reference to the $ k_{T} $ resolution.

Definition at line 242 of file Lepjets_Event.cc.

References _kt_res.

249 {
250  return _kt_res;
251 }
Lepjets_Event_Lep & hitfit::Lepjets_Event::lep ( std::vector< Lepjets_Event_Lep >::size_type  i)

Return a reference to lepton at index position i.

Parameters
iThe lepton index position.

Definition at line 142 of file Lepjets_Event.cc.

References _leps, assert(), and i.

Referenced by TtSemiLepHitFitProducer< LeptonCollection >::produce().

152 {
153  assert (i < _leps.size());
154  return _leps[i];
155 }
int i
Definition: DBlmapReader.cc:9
assert(m_qm.get())
std::vector< Lepjets_Event_Lep > _leps
const Lepjets_Event_Lep & hitfit::Lepjets_Event::lep ( std::vector< Lepjets_Event_Lep >::size_type  i) const

Return a constant reference to lepton at index position i.

Parameters
iThe lepton index position.

Definition at line 174 of file Lepjets_Event.cc.

References assert(), and i.

184 {
185  assert (i < _leps.size());
186  return _leps[i];
187 }
int i
Definition: DBlmapReader.cc:9
assert(m_qm.get())
std::vector< Lepjets_Event_Lep > _leps
Fourvec & hitfit::Lepjets_Event::met ( )
const Fourvec & hitfit::Lepjets_Event::met ( ) const

Return a constant reference to the missing transverse energy.

Definition at line 218 of file Lepjets_Event.cc.

References _met.

225 {
226  return _met;
227 }
std::vector< Lepjets_Event_Jet >::size_type hitfit::Lepjets_Event::njets ( ) const

Return the number of jets in the event.

Definition at line 130 of file Lepjets_Event.cc.

References _jets.

Referenced by hitfit::Top_Fit::fit(), and set_jet_types().

137 {
138  return _jets.size ();
139 }
std::vector< Lepjets_Event_Jet > _jets
std::vector< Lepjets_Event_Lep >::size_type hitfit::Lepjets_Event::nleps ( ) const

Return the number of leptons in the event.

Definition at line 118 of file Lepjets_Event.cc.

References _leps.

125 {
126  return _leps.size ();
127 }
std::vector< Lepjets_Event_Lep > _leps
const int & hitfit::Lepjets_Event::runnum ( ) const

Return a constant reference to the run number.

Definition at line 70 of file Lepjets_Event.cc.

References _runnum.

77 {
78  return _runnum;
79 }
int & hitfit::Lepjets_Event::runnum ( )

Return a reference to the run number.

Definition at line 82 of file Lepjets_Event.cc.

References _runnum.

89 {
90  return _runnum;
91 }
bool hitfit::Lepjets_Event::set_jet_types ( const std::vector< int > &  _jet_types)

Set the jet types in the event.

Definition at line 464 of file Lepjets_Event.cc.

References _jets, hitfit::hadw1_label, hitfit::hadw2_label, i, jet(), njets(), lumiQTWidget::t, and hitfit::Lepjets_Event_Lep::type().

Referenced by hitfit::RunHitFit< AElectron, AMuon, AJet, AMet >::FitAllPermutation().

469 {
470  if (_jets.size() != _jet_types.size()) {
471  return false;
472  }
473  bool saw_hadw1 = false;
475  int t = _jet_types[i];
476  if (t == hadw1_label) {
477  if (saw_hadw1)
478  t = hadw2_label;
479  saw_hadw1 = true;
480  }
481  jet (i).type() = t;
482  }
483  return true;
484 }
int i
Definition: DBlmapReader.cc:9
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.
int & type()
Return a reference to the type code.
std::vector< Lepjets_Event_Jet >::size_type njets() const
Return the number of jets in the event.
std::vector< Lepjets_Event_Jet > _jets
void hitfit::Lepjets_Event::setMC ( bool  isMC)

Set the Monte Carlo flag.

Definition at line 290 of file Lepjets_Event.cc.

References _isMC, and isMC().

297 {
298  _isMC = isMC;
299 }
bool isMC() const
Return the Monte Carlo flag.
void hitfit::Lepjets_Event::smear ( CLHEP::HepRandomEngine &  engine,
bool  smear_dir = false 
)

Smear the objects in the event according to their resolutions.

Parameters
engineThe underlying random number generator.
smear_dirIf TRUE, also smear the object's direction.
If FALSE, then only smear the magnitude of three-momentum.

Definition at line 411 of file Lepjets_Event.cc.

References _jets, _kt_res, _leps, _met, i, kt(), hitfit::Resolution::pick(), and X.

Referenced by hitfit::gentop(), and hitfit::gentth().

419 {
420  Fourvec before, after;
421  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
422  before += _leps[i].p();
423  _leps[i].smear (engine, smear_dir);
424  after += _leps[i].p();
425  }
426  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
427  before += _jets[i].p();
428  _jets[i].smear (engine, smear_dir);
429  after += _jets[i].p();
430  }
431 
432  Fourvec kt = _met + before;
434  kt(Fourvec::Y) = _kt_res.pick (kt(Fourvec::Y), kt(Fourvec::Y), engine);
435  _met = kt - after;
436 }
int i
Definition: DBlmapReader.cc:9
Fourvec kt() const
Return the sum of all objects&#39; four-momentum and missing transverse energy.
#define X(str)
Definition: MuonsGrabber.cc:48
std::vector< Lepjets_Event_Lep > _leps
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:227
uint16_t size_type
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:57
std::vector< Lepjets_Event_Jet > _jets
void hitfit::Lepjets_Event::sort ( )

Sort objects in the event according to their transverse momentum $ p_{T} $ .

Definition at line 439 of file Lepjets_Event.cc.

References _jets, and _leps.

443 {
444  std::stable_sort (_leps.begin(), _leps.end(), not2 (less<Lepjets_Event_Lep> ()));
445  std::stable_sort (_jets.begin(), _jets.end(), not2 (less<Lepjets_Event_Lep> ()));
446 }
std::vector< Lepjets_Event_Lep > _leps
std::vector< Lepjets_Event_Jet > _jets
Fourvec hitfit::Lepjets_Event::sum ( int  type) const

Return the sum of all objects' four-momentum which have a particular type.

Parameters
typeThe type code of the objects to be summed up.

Definition at line 349 of file Lepjets_Event.cc.

References i, dbtoconf::out, and AlCaHLTBitMon_ParallelJobs::p.

Referenced by hitfit::Top_Decaykin::hadt(), hitfit::Top_Decaykin::hadw(), hitfit::Top_Decaykin::hadw1(), hitfit::Top_Decaykin::hadw2(), hitfit::Top_Decaykin::lept(), and hitfit::Top_Decaykin::solve_nu_tmass().

359 {
360  Fourvec out;
362  if (_leps[i].type() == type)
363  out += _leps[i].p();
365  if (_jets[i].type() == type)
366  out += _jets[i].p();
367  return out;
368 }
type
Definition: HCALResponse.h:21
int i
Definition: DBlmapReader.cc:9
std::vector< Lepjets_Event_Lep > _leps
uint16_t size_type
CLHEP::HepLorentzVector Fourvec
Typedef for a HepLorentzVector.
Definition: fourvec.h:57
tuple out
Definition: dbtoconf.py:99
std::vector< Lepjets_Event_Jet > _jets
void hitfit::Lepjets_Event::trimjets ( std::vector< Lepjets_Event_Jet >::size_type  n)

Remove all but the first n jets.

Parameters
nThe number of jets to keep.

Definition at line 554 of file Lepjets_Event.cc.

References _jets, and gen::n.

561 {
562  if (n >= _jets.size())
563  return;
564  _jets.erase (_jets.begin() + n, _jets.end());
565 }
std::vector< Lepjets_Event_Jet > _jets
double hitfit::Lepjets_Event::zvertex ( ) const

Return the value of z-vertex.

Definition at line 254 of file Lepjets_Event.cc.

References _zvertex.

261 {
262  return _zvertex;
263 }
double & hitfit::Lepjets_Event::zvertex ( )

Return a reference to the value of z-vertex.

Definition at line 266 of file Lepjets_Event.cc.

References _zvertex.

273 {
274  return _zvertex;
275 }

Member Data Documentation

double hitfit::Lepjets_Event::_dlb
private

The low-bias (LB) discriminant.

Definition at line 383 of file Lepjets_Event.h.

Referenced by dlb().

double hitfit::Lepjets_Event::_dnn
private

The neural network (NN) discriminant.

Definition at line 388 of file Lepjets_Event.h.

Referenced by dnn().

int hitfit::Lepjets_Event::_evnum
private

The event number.

Definition at line 378 of file Lepjets_Event.h.

Referenced by evnum().

bool hitfit::Lepjets_Event::_isMC
private

The Monte Calro flag.

Definition at line 368 of file Lepjets_Event.h.

Referenced by isMC(), and setMC().

std::vector<Lepjets_Event_Jet> hitfit::Lepjets_Event::_jets
private

The list of jets in the event.

Definition at line 347 of file Lepjets_Event.h.

Referenced by add_jet(), cut_jets(), jet(), jet_permutation(), jet_types(), kt(), njets(), set_jet_types(), smear(), sort(), and trimjets().

Resolution hitfit::Lepjets_Event::_kt_res
private

The $ k_{T} $ resolution.

Definition at line 358 of file Lepjets_Event.h.

Referenced by kt_res(), and smear().

std::vector<Lepjets_Event_Lep> hitfit::Lepjets_Event::_leps
private

The list of leptons in the event.

Definition at line 342 of file Lepjets_Event.h.

Referenced by add_lep(), cut_leps(), kt(), lep(), nleps(), smear(), and sort().

Fourvec hitfit::Lepjets_Event::_met
private

Missing transverse energy.

Definition at line 353 of file Lepjets_Event.h.

Referenced by kt(), met(), and smear().

int hitfit::Lepjets_Event::_runnum
private

The run number.

Definition at line 373 of file Lepjets_Event.h.

Referenced by runnum().

double hitfit::Lepjets_Event::_zvertex
private

The $ z- $ vertex of the event.

Definition at line 363 of file Lepjets_Event.h.

Referenced by zvertex().