CMS 3D CMS Logo

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>

List of all members.

Public Member Functions

void add_jet (const Lepjets_Event_Jet &jet)
 Add a new jet to the event.
void add_lep (const Lepjets_Event_Lep &lep)
 Add a new lepton to the event.
int cut_jets (double pt_cut, double eta_cut)
 Remove jets which fail transverse momentum $ p_{T} $ and pseudorapidity $ \eta $ cut.
int cut_leps (double pt_cut, double eta_cut)
 Remove leptons which fail transverse momentum $ p_{T} $ and pseudorapidity $ \eta $ cut.
double dlb () const
 Return the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment).
double & dlb ()
 Return a reference to the value of low-bias (LB) discriminant (Irrelevant for non-D0 experiment).
double & dnn ()
 Return a reference to the value of neural network (NN) discriminant (Irrelevant for non-D0 experiment).
double dnn () const
 Return a the value of neural network (NN) discriminant (Irrelevant for non-D0 experiment).
std::ostream & dump (std::ostream &s, bool full=false) const
 Print the content of this object.
int & evnum ()
 Return a reference to the event number.
const int & evnum () const
 Return a constant reference to the event number.
bool isMC () const
 Return the Monte Carlo flag.
Lepjets_Event_Jetjet (std::vector< Lepjets_Event_Jet >::size_type i)
 Return a reference to jet at index position i.
const Lepjets_Event_Jetjet (std::vector< Lepjets_Event_Jet >::size_type i) const
 Return a constant reference to jet at index position i.
std::string 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.

std::vector< int > jet_types () const
 Return the jet types in the event.
Fourvec kt () const
 Return the sum of all objects' four-momentum and missing transverse energy.
Resolutionkt_res ()
 Return a reference to the $ k_{T} $ resolution.
const Resolutionkt_res () const
 Return a const reference to the $ k_{T} $ resolution.
Lepjets_Event_Leplep (std::vector< Lepjets_Event_Lep >::size_type i)
 Return a reference to lepton at index position i.
const Lepjets_Event_Leplep (std::vector< Lepjets_Event_Lep >::size_type i) const
 Return a constant reference to lepton at index position i.
 Lepjets_Event (int runnum, int evnum)
 Constructor.
const Fourvecmet () const
 Return a constant reference to the missing transverse energy.
Fourvecmet ()
 Return a reference to the missing transverse energy.
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.
const int & runnum () const
 Return a constant reference to the run number.
int & runnum ()
 Return a reference to the run number.
bool set_jet_types (const std::vector< int > &)
 Set the jet types in the event.
void setMC (bool isMC)
 Set the Monte Carlo flag.
void smear (CLHEP::HepRandomEngine &engine, bool smear_dir=false)
 Smear the objects in the event according to their resolutions.
void sort ()
 Sort objects in the event according to their transverse momentum $ p_{T} $ .
Fourvec sum (int type) const
 Return the sum of all objects' four-momentum which have a particular type.
void trimjets (std::vector< Lepjets_Event_Jet >::size_type n)
 Remove all but the first n jets.
double zvertex () const
 Return the value of z-vertex.
double & zvertex ()
 Return a reference to the value of z-vertex.

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

  : _zvertex (0),
    _isMC(false),
    _runnum (runnum),
    _evnum (evnum),
    _dlb (-1),
    _dnn (-1)
{
}

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

References _jets.

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

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

References _leps.

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

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

References _jets.

{
  _jets.erase (remove_if (_jets.begin(), _jets.end(),
                          Lepjets_Event_Cutter (pt_cut, eta_cut)),
               _jets.end ());
  return _jets.size ();
}
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 517 of file Lepjets_Event.cc.

References _leps.

{
  _leps.erase (remove_if (_leps.begin(), _leps.end(),
                          Lepjets_Event_Cutter (pt_cut, eta_cut)),
               _leps.end ());
  return _leps.size ();
}
double hitfit::Lepjets_Event::dlb ( ) const

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

Definition at line 302 of file Lepjets_Event.cc.

References _dlb.

{
  return _dlb;
}
double & hitfit::Lepjets_Event::dlb ( )

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

Definition at line 314 of file Lepjets_Event.cc.

References _dlb.

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

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

Definition at line 338 of file Lepjets_Event.cc.

References _dnn.

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

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

Definition at line 326 of file Lepjets_Event.cc.

References _dnn.

{
  return _dnn;
}
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 569 of file Lepjets_Event.cc.

References full, i, and alignCSCRings::s.

Referenced by hitfit::operator<<().

{
  s << "Run: " << _runnum << "  Event: " << _evnum << "\n";
  s << "Leptons:\n";
  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
    s << "  ";
    _leps[i].dump (s, full);
    s << "\n";
  }
  s << "Jets:\n";
  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
    s << "  ";
    _jets[i].dump (s, full);
    s << "\n";
  }
  s << "Missing Et: " << _met << "\n";
  if (_zvertex != 0)
    s << "z-vertex: " << _zvertex << "\n";
  return s;
}
int & hitfit::Lepjets_Event::evnum ( )

Return a reference to the event number.

Definition at line 107 of file Lepjets_Event.cc.

References _evnum.

{
  return _evnum;
}
const int & hitfit::Lepjets_Event::evnum ( ) const

Return a constant reference to the event number.

Definition at line 95 of file Lepjets_Event.cc.

References _evnum.

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

Return the Monte Carlo flag.

Definition at line 279 of file Lepjets_Event.cc.

References _isMC.

Referenced by setMC().

{
  return _isMC;
}
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 159 of file Lepjets_Event.cc.

References _jets, and i.

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

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

References i.

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

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

{
    std::string permutation;
    for (size_t jet = 0 ; jet != _jets.size() ; ++jet) {
        permutation = permutation + hitfit::jetTypeString(_jets[jet].type());
    }
    return permutation;
}
std::vector< int > hitfit::Lepjets_Event::jet_types ( ) const

Return the jet types in the event.

Definition at line 450 of file Lepjets_Event.cc.

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

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

{
  std::vector<int> ret;
  for (std::vector<Lepjets_Event_Jet>::size_type ijet =  0 ;
       ijet != _jets.size() ;
       ijet++) {
    ret.push_back(jet(ijet).type());
  }
  return ret;
}
Fourvec hitfit::Lepjets_Event::kt ( ) const

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

Definition at line 372 of file Lepjets_Event.cc.

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

Referenced by smear().

{
  Fourvec v = _met;
  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++)
    v += _leps[i].p();
  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++)
    v += _jets[i].p();
  return v;
}
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 243 of file Lepjets_Event.cc.

References _kt_res.

{
  return _kt_res;
}
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 143 of file Lepjets_Event.cc.

References _leps, and i.

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

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

References i.

{
  assert (i < _leps.size());
  return _leps[i];
}
Fourvec & hitfit::Lepjets_Event::met ( )
const Fourvec & hitfit::Lepjets_Event::met ( ) const

Return a constant reference to the missing transverse energy.

Definition at line 219 of file Lepjets_Event.cc.

References _met.

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

Return the number of jets in the event.

Definition at line 131 of file Lepjets_Event.cc.

References _jets.

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

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

Return the number of leptons in the event.

Definition at line 119 of file Lepjets_Event.cc.

References _leps.

{
  return _leps.size ();
}
int & hitfit::Lepjets_Event::runnum ( )

Return a reference to the run number.

Definition at line 83 of file Lepjets_Event.cc.

References _runnum.

{
  return _runnum;
}
const int & hitfit::Lepjets_Event::runnum ( ) const

Return a constant reference to the run number.

Definition at line 71 of file Lepjets_Event.cc.

References _runnum.

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

Set the jet types in the event.

Definition at line 465 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().

{
  if (_jets.size() != _jet_types.size()) {
    return false;
  }
  bool saw_hadw1 = false;
  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < njets(); i++) {
    int t = _jet_types[i];
    if (t == hadw1_label) {
      if (saw_hadw1)
        t = hadw2_label;
      saw_hadw1 = true;
    }
    jet (i).type() = t;
  }
  return true;
}
void hitfit::Lepjets_Event::setMC ( bool  isMC)

Set the Monte Carlo flag.

Definition at line 291 of file Lepjets_Event.cc.

References _isMC, and isMC().

{
  _isMC = isMC;
}
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 412 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().

{
  Fourvec before, after;
  for (std::vector<Lepjets_Event_Lep>::size_type i=0; i < _leps.size(); i++) {
    before += _leps[i].p();
    _leps[i].smear (engine, smear_dir);
    after += _leps[i].p();
  }
  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++) {
    before += _jets[i].p();
    _jets[i].smear (engine, smear_dir);
    after += _jets[i].p();
  }

  Fourvec kt = _met + before;
  kt(Fourvec::X) = _kt_res.pick (kt(Fourvec::X), kt(Fourvec::X), engine);
  kt(Fourvec::Y) = _kt_res.pick (kt(Fourvec::Y), kt(Fourvec::Y), engine);
  _met = kt - after;
}
void hitfit::Lepjets_Event::sort ( )

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

Definition at line 440 of file Lepjets_Event.cc.

References _jets, and _leps.

{
  std::stable_sort (_leps.begin(), _leps.end(), not2 (less<Lepjets_Event_Lep> ()));
  std::stable_sort (_jets.begin(), _jets.end(), not2 (less<Lepjets_Event_Lep> ()));
}
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 350 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().

{
  Fourvec out;
  for (std::vector<Lepjets_Event_Lep>::size_type  i=0; i < _leps.size(); i++)
    if (_leps[i].type() == type)
      out += _leps[i].p();
  for (std::vector<Lepjets_Event_Jet>::size_type i=0; i < _jets.size(); i++)
    if (_jets[i].type() == type)
      out += _jets[i].p();
  return out;
}
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 555 of file Lepjets_Event.cc.

References _jets, and n.

{
  if (n >= _jets.size())
    return;
  _jets.erase (_jets.begin() + n, _jets.end());
}
double & hitfit::Lepjets_Event::zvertex ( )

Return a reference to the value of z-vertex.

Definition at line 267 of file Lepjets_Event.cc.

References _zvertex.

{
  return _zvertex;
}
double hitfit::Lepjets_Event::zvertex ( ) const

Return the value of z-vertex.

Definition at line 255 of file Lepjets_Event.cc.

References _zvertex.

{
  return _zvertex;
}

Member Data Documentation

double hitfit::Lepjets_Event::_dlb [private]

The low-bias (LB) discriminant.

Definition at line 384 of file Lepjets_Event.h.

Referenced by dlb().

double hitfit::Lepjets_Event::_dnn [private]

The neural network (NN) discriminant.

Definition at line 389 of file Lepjets_Event.h.

Referenced by dnn().

The event number.

Definition at line 379 of file Lepjets_Event.h.

Referenced by evnum().

The Monte Calro flag.

Definition at line 369 of file Lepjets_Event.h.

Referenced by isMC(), and setMC().

The list of jets in the event.

Definition at line 348 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().

The $ k_{T} $ resolution.

Definition at line 359 of file Lepjets_Event.h.

Referenced by kt_res(), and smear().

The list of leptons in the event.

Definition at line 343 of file Lepjets_Event.h.

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

Missing transverse energy.

Definition at line 354 of file Lepjets_Event.h.

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

The run number.

Definition at line 374 of file Lepjets_Event.h.

Referenced by runnum().

The $ z- $ vertex of the event.

Definition at line 364 of file Lepjets_Event.h.

Referenced by zvertex().