CMS 3D CMS Logo

TtEvent.h
Go to the documentation of this file.
1 #ifndef TopObjects_TtEvent_h
2 #define TopObjects_TtEvent_h
3 
4 #include <vector>
5 #include <string>
6 
11 
23 class TtEvent {
24 public:
26  enum HypoClassKey {
36  };
38  typedef std::pair<reco::CompositeCandidate, std::vector<int> > HypoCombPair;
41  const char* label;
43  };
44 
45 protected:
48 
49 public:
51  TtEvent(){};
53  virtual ~TtEvent(){};
54 
56  std::pair<WDecay::LeptonType, WDecay::LeptonType> lepDecays() const { return lepDecays_; }
59  const reco::CompositeCandidate& eventHypo(const HypoClassKey& key, const unsigned& cmb = 0) const {
60  return (evtHyp_.find(key)->second)[cmb].first;
61  };
63  const edm::RefProd<TtGenEvent>& genEvent() const { return genEvt_; };
64 
68  bool isHypoClassAvailable(const HypoClassKey& key) const { return (evtHyp_.find(key) != evtHyp_.end()); };
69  // check if hypothesis 'cmb' is available within the hypothesis class
70  bool isHypoAvailable(const std::string& key, const unsigned& cmb = 0) const {
72  };
74  bool isHypoAvailable(const HypoClassKey& key, const unsigned& cmb = 0) const {
75  return isHypoClassAvailable(key) ? (cmb < evtHyp_.find(key)->second.size()) : false;
76  };
78  bool isHypoValid(const std::string& key, const unsigned& cmb = 0) const {
80  };
82  bool isHypoValid(const HypoClassKey& key, const unsigned& cmb = 0) const {
83  return isHypoAvailable(key, cmb) ? !eventHypo(key, cmb).roles().empty() : false;
84  };
86  unsigned int numberOfAvailableHypoClasses() const { return evtHyp_.size(); };
88  unsigned int numberOfAvailableHypos(const std::string& key) const {
90  };
92  unsigned int numberOfAvailableHypos(const HypoClassKey& key) const {
93  return isHypoAvailable(key) ? evtHyp_.find(key)->second.size() : 0;
94  };
98  };
101  return (isHypoAvailable(key) ? nJetsConsidered_.find(key)->second : -1);
102  };
104  std::vector<int> jetLeptonCombination(const std::string& key, const unsigned& cmb = 0) const {
106  };
108  std::vector<int> jetLeptonCombination(const HypoClassKey& key, const unsigned& cmb = 0) const {
109  return (evtHyp_.find(key)->second)[cmb].second;
110  };
112  double genMatchSumPt(const unsigned& cmb = 0) const {
113  return (cmb < genMatchSumPt_.size() ? genMatchSumPt_[cmb] : -1.);
114  };
116  double genMatchSumDR(const unsigned& cmb = 0) const {
117  return (cmb < genMatchSumDR_.size() ? genMatchSumDR_[cmb] : -1.);
118  };
120  std::string mvaMethod() const { return mvaMethod_; }
122  double mvaDisc(const unsigned& cmb = 0) const { return (cmb < mvaDisc_.size() ? mvaDisc_[cmb] : -1.); }
124  double fitChi2(const unsigned& cmb = 0) const { return (cmb < fitChi2_.size() ? fitChi2_[cmb] : -1.); }
126  double hitFitChi2(const unsigned& cmb = 0) const { return (cmb < hitFitChi2_.size() ? hitFitChi2_[cmb] : -1.); }
128  double fitProb(const unsigned& cmb = 0) const { return (cmb < fitProb_.size() ? fitProb_[cmb] : -1.); }
130  double hitFitProb(const unsigned& cmb = 0) const { return (cmb < hitFitProb_.size() ? hitFitProb_[cmb] : -1.); }
132  double hitFitMT(const unsigned& cmb = 0) const { return (cmb < hitFitMT_.size() ? hitFitMT_[cmb] : -1.); }
134  double hitFitSigMT(const unsigned& cmb = 0) const { return (cmb < hitFitSigMT_.size() ? hitFitSigMT_[cmb] : -1.); }
136  int correspondingHypo(const std::string& key1, const unsigned& hyp1, const std::string& key2) const {
138  };
140  int correspondingHypo(const HypoClassKey& key1, const unsigned& hyp1, const HypoClassKey& key2) const;
141 
143  const reco::Candidate* topPair(const std::string& key, const unsigned& cmb = 0) const {
144  return topPair(hypoClassKeyFromString(key), cmb);
145  };
147  const reco::Candidate* topPair(const HypoClassKey& key, const unsigned& cmb = 0) const {
148  return !isHypoValid(key, cmb) ? nullptr : (reco::Candidate*)&eventHypo(key, cmb);
149  };
151  const math::XYZTLorentzVector* topPair() const { return (!genEvt_ ? nullptr : this->genEvent()->topPair()); };
152 
154  void setLepDecays(const WDecay::LeptonType& lepDecTop1, const WDecay::LeptonType& lepDecTop2) {
155  lepDecays_ = std::make_pair(lepDecTop1, lepDecTop2);
156  };
160  void addEventHypo(const HypoClassKey& key, const HypoCombPair& hyp) { evtHyp_[key].push_back(hyp); };
162  void setNumberOfConsideredJets(const HypoClassKey& key, const unsigned int nJets) { nJetsConsidered_[key] = nJets; };
164  void setGenMatchSumPt(const std::vector<double>& val) { genMatchSumPt_ = val; };
166  void setGenMatchSumDR(const std::vector<double>& val) { genMatchSumDR_ = val; };
170  void setMvaDiscriminators(const std::vector<double>& val) { mvaDisc_ = val; };
172  void setFitChi2(const std::vector<double>& val) { fitChi2_ = val; };
174  void setHitFitChi2(const std::vector<double>& val) { hitFitChi2_ = val; };
176  void setFitProb(const std::vector<double>& val) { fitProb_ = val; };
178  void setHitFitProb(const std::vector<double>& val) { hitFitProb_ = val; };
180  void setHitFitMT(const std::vector<double>& val) { hitFitMT_ = val; };
182  void setHitFitSigMT(const std::vector<double>& val) { hitFitSigMT_ = val; };
183 
184 protected:
186  std::pair<WDecay::LeptonType, WDecay::LeptonType> lepDecays_;
191  std::map<HypoClassKey, std::vector<HypoCombPair> > evtHyp_;
193  std::map<HypoClassKey, int> nJetsConsidered_;
194 
196  std::vector<double> fitChi2_;
197  std::vector<double> hitFitChi2_;
199  std::vector<double> fitProb_;
200  std::vector<double> hitFitProb_;
202  std::vector<double> hitFitMT_;
203  std::vector<double> hitFitSigMT_;
205  std::vector<double> genMatchSumPt_;
207  std::vector<double> genMatchSumDR_;
211  std::vector<double> mvaDisc_;
212 };
213 
214 #endif
unsigned int numberOfAvailableHypoClasses() const
return number of available hypothesis classes
Definition: TtEvent.h:86
unsigned int numberOfAvailableHypos(const std::string &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:88
void setMvaDiscriminators(const std::vector< double > &val)
set mva discriminant values of kMVADisc hypothesis
Definition: TtEvent.h:170
std::map< HypoClassKey, int > nJetsConsidered_
number of jets considered when building the hypotheses
Definition: TtEvent.h:193
HypoClassKey
supported classes of event hypotheses
Definition: TtEvent.h:26
void setHitFitSigMT(const std::vector< double > &val)
set fitted top mass uncertainty of kHitFit hypothesis
Definition: TtEvent.h:182
std::vector< double > genMatchSumDR_
result of gen match
Definition: TtEvent.h:207
TtEvent()
empty constructor
Definition: TtEvent.h:51
double hitFitProb(const unsigned &cmb=0) const
return the hitfit probability of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:130
const reco::Candidate * topPair(const HypoClassKey &key, const unsigned &cmb=0) const
get combined 4-vector of top and topBar of the given hypothesis
Definition: TtEvent.h:147
double fitChi2(const unsigned &cmb=0) const
return the chi2 of the kinematic fit of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:124
void setMvaMethod(const std::string &name)
set label of mva method for kMVADisc hypothesis
Definition: TtEvent.h:168
std::string mvaMethod() const
return the label of the mva method in use for the jet parton association (if kMVADisc is not availabl...
Definition: TtEvent.h:120
std::vector< double > hitFitChi2_
Definition: TtEvent.h:197
void setGenMatchSumDR(const std::vector< double > &val)
set sum dr of kGenMatch hypothesis
Definition: TtEvent.h:166
void setFitChi2(const std::vector< double > &val)
set chi2 of kKinFit hypothesis
Definition: TtEvent.h:172
const math::XYZTLorentzVector * topPair() const
get combined 4-vector of top and topBar from the TtGenEvent
Definition: TtEvent.h:151
void setLepDecays(const WDecay::LeptonType &lepDecTop1, const WDecay::LeptonType &lepDecTop2)
set leptonic decay channels
Definition: TtEvent.h:154
int numberOfConsideredJets(const HypoClassKey &key) const
return number of jets that were considered when building a given hypothesis
Definition: TtEvent.h:100
std::pair< reco::CompositeCandidate, std::vector< int > > HypoCombPair
pair of hypothesis and lepton jet combinatorics for a given hypothesis
Definition: TtEvent.h:38
std::vector< double > mvaDisc_
MVA discriminants.
Definition: TtEvent.h:211
const edm::RefProd< TtGenEvent > & genEvent() const
get TtGenEvent
Definition: TtEvent.h:63
int correspondingHypo(const std::string &key1, const unsigned &hyp1, const std::string &key2) const
return the hypothesis in hypothesis class &#39;key2&#39;, which corresponds to hypothesis &#39;hyp1&#39; in hypothesi...
Definition: TtEvent.h:136
std::pair< WDecay::LeptonType, WDecay::LeptonType > lepDecays_
leptonic decay channels
Definition: TtEvent.h:182
U second(std::pair< T, U > const &p)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
char const * label
virtual ~TtEvent()
default destructor
Definition: TtEvent.h:53
bool isHypoAvailable(const std::string &key, const unsigned &cmb=0) const
Definition: TtEvent.h:70
std::vector< double > fitProb_
result of kinematic fit
Definition: TtEvent.h:199
std::string mvaMethod_
label of the MVA method
Definition: TtEvent.h:209
void setGenMatchSumPt(const std::vector< double > &val)
set sum pt of kGenMatch hypothesis
Definition: TtEvent.h:164
role_collection const & roles() const
get the roles
void setHitFitMT(const std::vector< double > &val)
set fitted top mass of kHitFit hypothesis
Definition: TtEvent.h:180
void setHitFitChi2(const std::vector< double > &val)
set chi2 of kHitFit hypothesis
Definition: TtEvent.h:174
std::vector< double > genMatchSumPt_
result of gen match
Definition: TtEvent.h:205
double hitFitChi2(const unsigned &cmb=0) const
return the hitfit chi2 of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:126
double genMatchSumDR(const unsigned &cmb=0) const
return the sum dr of the generator match if available; -1 else
Definition: TtEvent.h:116
bool isHypoClassAvailable(const HypoClassKey &key) const
check if hypothesis class &#39;key&#39; was added to the event structure
Definition: TtEvent.h:68
double hitFitMT(const unsigned &cmb=0) const
return the hitfit top mass of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:132
std::vector< double > hitFitSigMT_
Definition: TtEvent.h:203
void setFitProb(const std::vector< double > &val)
set fit probability of kKinFit hypothesis
Definition: TtEvent.h:176
void setGenEvent(const edm::Handle< TtGenEvent > &evt)
set TtGenEvent
Definition: TtEvent.h:158
HypoClassKey hypoClassKeyFromString(const std::string &label) const
return the corresponding enum value from a string
Definition: TtEvent.cc:15
std::vector< double > hitFitMT_
result of hitfit
Definition: TtEvent.h:202
a lightweight map for selection type string label and enum value
Definition: TtEvent.h:40
void addEventHypo(const HypoClassKey &key, const HypoCombPair &hyp)
add new hypotheses
Definition: TtEvent.h:160
bool isHypoValid(const HypoClassKey &key, const unsigned &cmb=0) const
check if hypothesis &#39;cmb&#39; within the hypothesis class was valid; if not it lead to an empty Composite...
Definition: TtEvent.h:82
std::vector< double > hitFitProb_
Definition: TtEvent.h:200
bool isHypoClassAvailable(const std::string &key) const
check if hypothesis class &#39;key&#39; was added to the event structure
Definition: TtEvent.h:66
const reco::CompositeCandidate & eventHypo(const HypoClassKey &key, const unsigned &cmb=0) const
Definition: TtEvent.h:59
Base class to hold information for ttbar event interpretation.
Definition: TtEvent.h:23
unsigned int numberOfAvailableHypos(const HypoClassKey &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:92
edm::RefProd< TtGenEvent > genEvt_
reference to TtGenEvent (has to be kept in the event!)
Definition: TtEvent.h:188
void setNumberOfConsideredJets(const HypoClassKey &key, const unsigned int nJets)
set number of jets considered when building a given hypothesis
Definition: TtEvent.h:162
double fitProb(const unsigned &cmb=0) const
return the fit probability of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:128
std::vector< int > jetLeptonCombination(const std::string &key, const unsigned &cmb=0) const
return the vector of jet lepton combinatorics for a given hypothesis and class
Definition: TtEvent.h:104
bool isHypoValid(const std::string &key, const unsigned &cmb=0) const
check if hypothesis &#39;cmb&#39; within the hypothesis class was valid; if not it lead to an empty Composite...
Definition: TtEvent.h:78
std::vector< double > fitChi2_
result of kinematic fit
Definition: TtEvent.h:196
int numberOfConsideredJets(const std::string &key) const
return number of jets that were considered when building a given hypothesis
Definition: TtEvent.h:96
std::map< HypoClassKey, std::vector< HypoCombPair > > evtHyp_
Definition: TtEvent.h:191
double genMatchSumPt(const unsigned &cmb=0) const
return the sum pt of the generator match if available; -1 else
Definition: TtEvent.h:112
const reco::Candidate * topPair(const std::string &key, const unsigned &cmb=0) const
get combined 4-vector of top and topBar of the given hypothesis
Definition: TtEvent.h:143
std::vector< int > jetLeptonCombination(const HypoClassKey &key, const unsigned &cmb=0) const
return the vector of jet lepton combinatorics for a given hypothesis and class
Definition: TtEvent.h:108
bool isHypoAvailable(const HypoClassKey &key, const unsigned &cmb=0) const
check if hypothesis &#39;cmb&#39; is available within the hypothesis class
Definition: TtEvent.h:74
double mvaDisc(const unsigned &cmb=0) const
return the mva discriminant value of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:122
void setHitFitProb(const std::vector< double > &val)
set fit probability of kHitFit hypothesis
Definition: TtEvent.h:178
double hitFitSigMT(const unsigned &cmb=0) const
return the hitfit top mass uncertainty of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:134
std::pair< WDecay::LeptonType, WDecay::LeptonType > lepDecays() const
get leptonic decay channels
Definition: TtEvent.h:56