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};
24  enum Operator {kAdd, kMult};
25 }
26 
38 
39  public:
40 
44  TtSemiLepJetComb(const std::vector<pat::Jet>&, const std::vector<int>&, const math::XYZTLorentzVector&, const pat::MET&);
47 
49  double topVar(JetComb::DecayType decay, JetComb::VarType var) const;
51  double wBosonVar(JetComb::DecayType decay, JetComb::VarType var) const;
53  double bQuarkVar(JetComb::DecayType decay, JetComb::VarType var) const;
55  double lightQVar(bool qbar, JetComb::VarType var) const;
57  double leptonVar(JetComb::VarType var) const;
59  double neutrinoVar(JetComb::VarType var) const;
60 
62  double compareHadTopLepTop(JetComb::CompType comp) const;
64  double compareHadWLepW(JetComb::CompType comp) const;
66  double compareHadBLepB(JetComb::CompType comp) const;
68  double compareLightQuarks(JetComb::CompType comp) const;
70  double compareLeptonNeutrino(JetComb::CompType comp) const;
72  double compareTopW(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const;
74  double compareTopB(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const;
76  double compareWB(JetComb::DecayType dec1, JetComb::DecayType dec2, JetComb::CompType comp) const;
78  double compareTopLepton(JetComb::DecayType decay, JetComb::CompType comp) const;
80  double compareTopNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const;
82  double compareWLepton(JetComb::DecayType decay, JetComb::CompType comp) const;
84  double compareWNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const;
86  double compareBLepton(JetComb::DecayType decay, JetComb::CompType comp) const;
88  double compareBNeutrino(JetComb::DecayType decay, JetComb::CompType comp) const;
89 
93  double relativePtHadronicTop() const;
96  double bOverLightQPt() const;
97 
99  double bTag(JetComb::DecayType decay, JetComb::BTagAlgo algo) const { return bTag(bQuark(decay), algo); };
101  double combinedBTags(JetComb::BTagAlgo algo, JetComb::Operator op) const;
103  double combinedBTagsForLightQuarks(JetComb::BTagAlgo algo, JetComb::Operator op) const;
104 
106  double addUserVar(std::string key, double value) { return userVariables_[key]=value;};
108  double userVar(const std::string& key) const {
109  return (userVariables_.find(key)!=userVariables_.end() ? userVariables_.find(key)->second : -9999.);};
110 
111  private:
112 
114  void deduceMothers();
116  double bTag(const pat::Jet& jet, JetComb::BTagAlgo algo) const;
118  double lightQVar(JetComb::VarType var) const { return lightQVar(false, var); };
120  const pat::Jet& lightQ(bool qbar=false) const { return (qbar ? hadQBarJet_ : hadQJet_); }
122  const pat::Jet& bQuark(JetComb::DecayType decay) const { return (decay==JetComb::kHad ? hadBJet_ : lepBJet_); }
124  const math::XYZTLorentzVector& wBoson(JetComb::DecayType decay) const { return (decay==JetComb::kHad ? hadW_ : lepW_); }
126  const math::XYZTLorentzVector& top(JetComb::DecayType decay) const { return (decay==JetComb::kHad ? hadTop_ : lepTop_); }
127 
128  private:
129 
151  std::map<std::string, double> userVariables_;
152 };
153 
154 #endif
Analysis-level MET class.
Definition: MET.h:43
math::XYZTLorentzVector lepton_
lepton 4-vector
math::XYZTLorentzVector hadTop_
hadronic top 4-vector
CompType
supported comparison types
math::XYZTLorentzVector lepW_
leptonic W 4-vector
math::XYZTLorentzVector lepTop_
leptonic top 4-vector
std::map< std::string, double > userVariables_
map for user defined variables
const math::XYZTLorentzVector & top(JetComb::DecayType decay) const
return leptonic or hadronic top candidate depending on argument
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Operator
operators for combining variables
const pat::Jet & lightQ(bool qbar=false) const
return lightQ or lightQBar candidate depending on argument
pat::Jet hadBJet_
hadronic b jet
Definition: value.py:1
double lightQVar(JetComb::VarType var) const
light quark candidate variable with default on q and not qbar
pat::Jet hadQJet_
lightQ jet
pat::Jet lepBJet_
leptonic b jet
double addUserVar(std::string key, double value)
add an arbitary user defined variable with given key and value
math::XYZTLorentzVector hadW_
hadronic W 4-vector
double userVar(const std::string &key) const
receive user defined variable value with given key
Analysis-level calorimeter jet class.
Definition: Jet.h:78
VarType
supported std single variable types
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 & wBoson(JetComb::DecayType decay) const
return leptonic or hadronic W candidate depending on argument
const pat::Jet & bQuark(JetComb::DecayType decay) const
return leptonic or hadronic b candidate depending on argument
double bTag(JetComb::DecayType decay, JetComb::BTagAlgo algo) const
b-tag value of a b candidate
pat::MET neutrino_
neutrino candidate