CMS 3D CMS Logo

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 50 of file Lepjets_Event.cc.

58  : _zvertex (0),
59  _isMC(false),
60  _runnum (runnum),
61  _evnum (evnum),
62  _dlb (-1),
63  _dnn (-1)
64 {
65 }
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 397 of file Lepjets_Event.cc.

References _jets.

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

404 {
405  _jets.push_back (jet);
406 }
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 385 of file Lepjets_Event.cc.

References _leps.

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

392 {
393  _leps.push_back (lep);
394 }
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 533 of file Lepjets_Event.cc.

References _jets.

544 {
545  _jets.erase (remove_if (_jets.begin(), _jets.end(),
546  Lepjets_Event_Cutter (pt_cut, eta_cut)),
547  _jets.end ());
548  return _jets.size ();
549 }
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 514 of file Lepjets_Event.cc.

References _leps.

525 {
526  _leps.erase (remove_if (_leps.begin(), _leps.end(),
527  Lepjets_Event_Cutter (pt_cut, eta_cut)),
528  _leps.end ());
529  return _leps.size ();
530 }
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 311 of file Lepjets_Event.cc.

References _dlb.

318 {
319  return _dlb;
320 }
double hitfit::Lepjets_Event::dlb ( ) const

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

Definition at line 299 of file Lepjets_Event.cc.

References _dlb.

306 {
307  return _dlb;
308 }
double & hitfit::Lepjets_Event::dnn ( )

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

Definition at line 335 of file Lepjets_Event.cc.

References _dnn.

342 {
343  return _dnn;
344 }
double hitfit::Lepjets_Event::dnn ( ) const

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

Definition at line 323 of file Lepjets_Event.cc.

References _dnn.

330 {
331  return _dnn;
332 }
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 566 of file Lepjets_Event.cc.

References _evnum, _jets, _leps, _met, _runnum, _zvertex, full, mps_fire::i, and alignCSCRings::s.

Referenced by hitfit::operator<<().

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

Return a constant reference to the event number.

Definition at line 92 of file Lepjets_Event.cc.

References _evnum.

99 {
100  return _evnum;
101 }
int & hitfit::Lepjets_Event::evnum ( )

Return a reference to the event number.

Definition at line 104 of file Lepjets_Event.cc.

References _evnum.

111 {
112  return _evnum;
113 }
bool hitfit::Lepjets_Event::isMC ( ) const

Return the Monte Carlo flag.

Definition at line 276 of file Lepjets_Event.cc.

References _isMC.

Referenced by setMC().

283 {
284  return _isMC;
285 }
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 156 of file Lepjets_Event.cc.

References _jets, and mps_fire::i.

Referenced by hitfit::Constrained_Top::Constrained_Top(), hitfit::Constrained_Z::Constrained_Z(), hitfit::Top_Fit_Args::constrainer_args(), jet_permutation(), jet_types(), hitfit::Gentop_Args::kt_res_str(), TtSemiLepHitFitProducer< LeptonCollection >::produce(), and set_jet_types().

166 {
167  assert (i < _jets.size());
168  return _jets[i];
169 }
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 188 of file Lepjets_Event.cc.

References _jets, and mps_fire::i.

198 {
199  assert (i < _jets.size());
200  return _jets[i];
201 }
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 598 of file Lepjets_Event.cc.

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

607 {
608  std::string permutation;
609  for (size_t jet = 0 ; jet != _jets.size() ; ++jet) {
610  permutation = permutation + hitfit::jetTypeString(_jets[jet].type());
611  }
612  return permutation;
613 }
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 447 of file Lepjets_Event.cc.

References _jets, and jet().

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

451 {
452  std::vector<int> ret;
454  ijet != _jets.size() ;
455  ijet++) {
456  ret.push_back(jet(ijet).type());
457  }
458  return ret;
459 }
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 369 of file Lepjets_Event.cc.

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

Referenced by hitfit::Constrained_Top::Constrained_Top(), hitfit::Constrained_Z::Constrained_Z(), and smear().

375 {
376  Fourvec v = _met;
378  v += _leps[i].p();
380  v += _jets[i].p();
381  return v;
382 }
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 ( )
const Resolution & hitfit::Lepjets_Event::kt_res ( ) const

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

Definition at line 240 of file Lepjets_Event.cc.

References _kt_res.

247 {
248  return _kt_res;
249 }
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 140 of file Lepjets_Event.cc.

References _leps, and mps_fire::i.

Referenced by hitfit::Constrained_Top::Constrained_Top(), hitfit::Constrained_Z::Constrained_Z(), and TtSemiLepHitFitProducer< LeptonCollection >::produce().

150 {
151  assert (i < _leps.size());
152  return _leps[i];
153 }
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 172 of file Lepjets_Event.cc.

References _leps, and mps_fire::i.

182 {
183  assert (i < _leps.size());
184  return _leps[i];
185 }
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 216 of file Lepjets_Event.cc.

References _met.

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

Return the number of jets in the event.

Definition at line 128 of file Lepjets_Event.cc.

References _jets.

Referenced by hitfit::Constrained_Top::Constrained_Top(), hitfit::Constrained_Z::Constrained_Z(), hitfit::Top_Fit_Args::constrainer_args(), hitfit::Top_Fit::fit(), hitfit::Gentop_Args::kt_res_str(), and set_jet_types().

135 {
136  return _jets.size ();
137 }
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 116 of file Lepjets_Event.cc.

References _leps.

Referenced by hitfit::Constrained_Top::Constrained_Top(), and hitfit::Constrained_Z::Constrained_Z().

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

Return a constant reference to the run number.

Definition at line 68 of file Lepjets_Event.cc.

References _runnum.

75 {
76  return _runnum;
77 }
int & hitfit::Lepjets_Event::runnum ( )

Return a reference to the run number.

Definition at line 80 of file Lepjets_Event.cc.

References _runnum.

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

Set the jet types in the event.

Definition at line 462 of file Lepjets_Event.cc.

References _eta_cut, _jets, _pt_cut, hitfit::hadw1_label, hitfit::hadw2_label, mps_fire::i, jet(), checklumidiff::l, njets(), lumiQTWidget::t, and hitfit::Lepjets_Event_Lep::type().

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

467 {
468  if (_jets.size() != _jet_types.size()) {
469  return false;
470  }
471  bool saw_hadw1 = false;
473  int t = _jet_types[i];
474  if (t == hadw1_label) {
475  if (saw_hadw1)
476  t = hadw2_label;
477  saw_hadw1 = true;
478  }
479  jet (i).type() = t;
480  }
481  return true;
482 }
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 288 of file Lepjets_Event.cc.

References _isMC, and isMC().

295 {
296  _isMC = isMC;
297 }
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 409 of file Lepjets_Event.cc.

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

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

417 {
418  Fourvec before, after;
419  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
420  before += _leps[i].p();
421  _leps[i].smear (engine, smear_dir);
422  after += _leps[i].p();
423  }
424  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
425  before += _jets[i].p();
426  _jets[i].smear (engine, smear_dir);
427  after += _jets[i].p();
428  }
429 
430  Fourvec kt = _met + before;
433  _met = kt - after;
434 }
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 437 of file Lepjets_Event.cc.

References _jets, and _leps.

441 {
442  std::stable_sort (_leps.begin(), _leps.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
443  std::stable_sort (_jets.begin(), _jets.end(), std::not_fn(std::less<Lepjets_Event_Lep>()));
444 }
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 347 of file Lepjets_Event.cc.

References _jets, _leps, mps_fire::i, MillePedeFileConverter_cfg::out, and AlCaHLTBitMon_ParallelJobs::p.

Referenced by hitfit::Top_Fit_Args::constrainer_args(), 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().

357 {
358  Fourvec out;
360  if (_leps[i].type() == type)
361  out += _leps[i].p();
363  if (_jets[i].type() == type)
364  out += _jets[i].p();
365  return out;
366 }
type
Definition: HCALResponse.h:21
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
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 552 of file Lepjets_Event.cc.

References _jets, and gen::n.

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

Return the value of z-vertex.

Definition at line 252 of file Lepjets_Event.cc.

References _zvertex.

259 {
260  return _zvertex;
261 }
double & hitfit::Lepjets_Event::zvertex ( )

Return a reference to the value of z-vertex.

Definition at line 264 of file Lepjets_Event.cc.

References _zvertex.

271 {
272  return _zvertex;
273 }

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 dump(), and 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(), dump(), jet(), jet_permutation(), jet_types(), kt(), njets(), set_jet_types(), smear(), sort(), sum(), 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(), dump(), kt(), lep(), nleps(), smear(), sort(), and sum().

Fourvec hitfit::Lepjets_Event::_met
private

Missing transverse energy.

Definition at line 353 of file Lepjets_Event.h.

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

int hitfit::Lepjets_Event::_runnum
private

The run number.

Definition at line 373 of file Lepjets_Event.h.

Referenced by dump(), and runnum().

double hitfit::Lepjets_Event::_zvertex
private

The $ z- $ vertex of the event.

Definition at line 363 of file Lepjets_Event.h.

Referenced by dump(), and zvertex().