CMS 3D CMS Logo

TtSemiLepJetComb.h
Go to the documentation of this file.
1 #ifndef TtSemiLepJetComb_h
2 #define TtSemiLepJetComb_h
3 
4 #include <vector>
5 #include <string>
6 
7 #include "TMath.h"
8 
11 
12 namespace JetComb {
15  enum DecayType { kHad, kLep };
17  enum VarType { kMass, kPt, kEta, kPhi, kTheta };
21  enum BTagAlgo {
33  };
35  enum Operator { kAdd, kMult };
36 } // namespace JetComb
37 
49 public:
53  TtSemiLepJetComb(const std::vector<pat::Jet>&,
54  const std::vector<int>&,
56  const pat::MET&);
59 
67  double lightQVar(bool qbar, JetComb::VarType var) const;
69  double leptonVar(JetComb::VarType var) const;
71  double neutrinoVar(JetComb::VarType var) const;
72 
101 
105  double relativePtHadronicTop() const;
108  double bOverLightQPt() const;
109 
116 
118  double addUserVar(std::string key, double value) { return userVariables_[key] = value; };
120  double userVar(const std::string& key) const {
121  return (userVariables_.find(key) != userVariables_.end() ? userVariables_.find(key)->second : -9999.);
122  };
123 
124 private:
126  void deduceMothers();
128  double bTag(const pat::Jet& jet, JetComb::BTagAlgo algo) const;
130  double lightQVar(JetComb::VarType var) const { return lightQVar(false, var); };
132  const pat::Jet& lightQ(bool qbar = false) const { return (qbar ? hadQBarJet_ : hadQJet_); }
137  return (decay == JetComb::kHad ? hadW_ : lepW_);
138  }
141  return (decay == JetComb::kHad ? hadTop_ : lepTop_);
142  }
143 
144 private:
166  std::map<std::string, double> userVariables_;
167 };
168 
169 #endif
Analysis-level MET class.
Definition: MET.h:40
math::XYZTLorentzVector lepton_
lepton 4-vector
math::XYZTLorentzVector hadTop_
hadronic top 4-vector
double topVar(JetComb::DecayType decay, JetComb::VarType var) const
top quark candidate variable
CompType
supported comparison types
double compareTopW(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const
comparison between the top and the W candidate
math::XYZTLorentzVector lepW_
leptonic W 4-vector
double compareWB(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const
comparison between the W and the b candidate
math::XYZTLorentzVector lepTop_
leptonic top 4-vector
double leptonVar(JetComb::VarType var) const
lepton candidate variable
std::map< std::string, double > userVariables_
map for user defined variables
double wBosonVar(JetComb::DecayType decay, JetComb::VarType var) const
W boson candidate variable.
double compareTopNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the top and the neutrino candidate
double lightQVar(bool qbar, JetComb::VarType var) const
light quark candidate variable
~TtSemiLepJetComb()
default destructor
TtSemiLepJetComb()
emtpy constructor
double lightQVar(JetComb::VarType var) const
light quark candidate variable with default on q and not qbar
double neutrinoVar(JetComb::VarType var) const
neutrino candidate variable
double relativePtHadronicTop() const
double compareWLepton(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the W and the lepton candidate
const math::XYZTLorentzVector & wBoson(JetComb::DecayType decay) const
return leptonic or hadronic W candidate depending on argument
double compareTopB(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const
comparison between the top and the b candidate
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double compareHadWLepW(JetComb::CompType comp) const
comparison between the two W candidates
double compareHadTopLepTop(JetComb::CompType comp) const
comparison between the two top candidates
double userVar(const std::string &key) const
receive user defined variable value with given key
double combinedBTagsForLightQuarks(JetComb::BTagAlgo algo, JetComb::Operator op) const
combined b-tag values of the two light quark candidates
double bQuarkVar(JetComb::DecayType decay, JetComb::VarType var) const
b quark candidate variable
Operator
operators for combining variables
double compareLeptonNeutrino(JetComb::CompType comp) const
comparison between the lepton and the neutrino candidate
double compareLightQuarks(JetComb::CompType comp) const
comparison between the lightQ and the lightQBar candidate
pat::Jet hadBJet_
hadronic b jet
Definition: value.py:1
double compareBLepton(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the b and the lepton candidate
pat::Jet hadQJet_
lightQ jet
const pat::Jet & lightQ(bool qbar=false) const
return lightQ or lightQBar candidate depending on argument
pat::Jet lepBJet_
leptonic b jet
double addUserVar(std::string key, double value)
add an arbitary user defined variable with given key and value
double compareTopLepton(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the top and the lepton candidate
math::XYZTLorentzVector hadW_
hadronic W 4-vector
Analysis-level calorimeter jet class.
Definition: Jet.h:77
const pat::Jet & bQuark(JetComb::DecayType decay) const
return leptonic or hadronic b candidate depending on argument
VarType
supported std single variable types
double combinedBTags(JetComb::BTagAlgo algo, JetComb::Operator op) const
combined b-tag values of the two b candidates
pat::Jet hadQBarJet_
lightQBar jet
BTagAlgo
b-tagging algorithms
Common calculator class to keep multivariate analysis variables for jet combinations in semi-leptonic...
const math::XYZTLorentzVector & top(JetComb::DecayType decay) const
return leptonic or hadronic top candidate depending on argument
void deduceMothers()
reconstruct mother candidates from final state candidates
double compareBNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the b and the neutrino candidate
double bOverLightQPt() const
double compareWNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const
comparison between the W and the neutrino candidate
pat::MET neutrino_
neutrino candidate
double compareHadBLepB(JetComb::CompType comp) const
comparison between the two b candidates
double bTag(JetComb::DecayType decay, JetComb::BTagAlgo algo) const
b-tag value of a b candidate