CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 namespace edm{
24  class LogInfo;
25 }
26 
27 class TtEvent {
28 
29  public:
33  typedef std::pair<reco::CompositeCandidate, std::vector<int> > HypoCombPair;
34 
35  protected:
39  HypoClassKey hypoClassKeyFromString(const std::string& label) const;
40 
41  public:
43  TtEvent(){};
45  virtual ~TtEvent(){};
46 
48  std::pair<WDecay::LeptonType, WDecay::LeptonType> lepDecays() const { return lepDecays_; }
51  const reco::CompositeCandidate& eventHypo(const HypoClassKey& key, const unsigned& cmb=0) const { return (evtHyp_.find(key)->second)[cmb].first; };
53  const edm::RefProd<TtGenEvent>& genEvent() const { return genEvt_; };
54 
56  bool isHypoClassAvailable(const std::string& key) const { return isHypoClassAvailable( hypoClassKeyFromString(key) ); };
58  bool isHypoClassAvailable(const HypoClassKey& key) const { return (evtHyp_.find(key)!=evtHyp_.end()); };
59  // check if hypothesis 'cmb' is available within the hypothesis class
60  bool isHypoAvailable(const std::string& key, const unsigned& cmb=0) const { return isHypoAvailable( hypoClassKeyFromString(key), cmb ); };
62  bool isHypoAvailable(const HypoClassKey& key, const unsigned& cmb=0) const { return isHypoClassAvailable(key) ? (cmb<evtHyp_.find(key)->second.size()) : false; };
64  bool isHypoValid(const std::string& key, const unsigned& cmb=0) const { return isHypoValid( hypoClassKeyFromString(key), cmb ); };
66  bool isHypoValid(const HypoClassKey& key, const unsigned& cmb=0) const { return isHypoAvailable(key, cmb) ? !eventHypo(key, cmb).roles().empty() : false; };
68  unsigned int numberOfAvailableHypoClasses() const { return evtHyp_.size(); };
70  unsigned int numberOfAvailableHypos(const std::string& key) const { return numberOfAvailableHypos( hypoClassKeyFromString(key) ); };
72  unsigned int numberOfAvailableHypos(const HypoClassKey& key) const { return isHypoAvailable(key) ? evtHyp_.find(key)->second.size() : 0; };
74  int numberOfConsideredJets(const std::string& key) const { return numberOfConsideredJets(hypoClassKeyFromString(key) ); };
76  int numberOfConsideredJets(const HypoClassKey& key) const { return (isHypoAvailable(key) ? nJetsConsidered_.find(key)->second : -1); };
78  std::vector<int> jetLeptonCombination(const std::string& key, const unsigned& cmb=0) const { return jetLeptonCombination(hypoClassKeyFromString(key), cmb); };
80  std::vector<int> jetLeptonCombination(const HypoClassKey& key, const unsigned& cmb=0) const { return (evtHyp_.find(key)->second)[cmb].second; };
82  double genMatchSumPt(const unsigned& cmb=0) const { return (cmb<genMatchSumPt_.size() ? genMatchSumPt_[cmb] : -1.); };
84  double genMatchSumDR(const unsigned& cmb=0) const { return (cmb<genMatchSumDR_.size() ? genMatchSumDR_[cmb] : -1.); };
86  std::string mvaMethod() const { return mvaMethod_; }
88  double mvaDisc(const unsigned& cmb=0) const { return (cmb<mvaDisc_.size() ? mvaDisc_[cmb] : -1.); }
90  double fitChi2(const unsigned& cmb=0) const { return (cmb<fitChi2_.size() ? fitChi2_[cmb] : -1.); }
92  double hitFitChi2(const unsigned& cmb=0) const { return (cmb<hitFitChi2_.size() ? hitFitChi2_[cmb] : -1.); }
94  double fitProb(const unsigned& cmb=0) const { return (cmb<fitProb_.size() ? fitProb_[cmb] : -1.); }
96  double hitFitProb(const unsigned& cmb=0) const { return (cmb<hitFitProb_.size() ? hitFitProb_[cmb] : -1.); }
98  double hitFitMT(const unsigned& cmb=0) const { return (cmb<hitFitMT_.size() ? hitFitMT_[cmb] : -1.); }
100  double hitFitSigMT(const unsigned& cmb=0) const { return (cmb<hitFitSigMT_.size() ? hitFitSigMT_[cmb] : -1.); }
102  int correspondingHypo(const std::string& key1, const unsigned& hyp1, const std::string& key2) const { return correspondingHypo(hypoClassKeyFromString(key1), hyp1, hypoClassKeyFromString(key2) ); };
104  int correspondingHypo(const HypoClassKey& key1, const unsigned& hyp1, const HypoClassKey& key2) const;
105 
107  const reco::Candidate* topPair(const std::string& key, const unsigned& cmb=0) const { return topPair(hypoClassKeyFromString(key), cmb); };
109  const reco::Candidate* topPair(const HypoClassKey& key, const unsigned& cmb=0) const { return !isHypoValid(key,cmb) ? 0 : (reco::Candidate*)&eventHypo(key,cmb); };
111  const math::XYZTLorentzVector* topPair() const { return (!genEvt_ ? 0 : this->genEvent()->topPair()); };
112 
114  void printParticle(edm::LogInfo &log, const char* name, const reco::Candidate* cand) const;
115 
117  void setLepDecays(const WDecay::LeptonType& lepDecTop1, const WDecay::LeptonType& lepDecTop2) { lepDecays_=std::make_pair(lepDecTop1, lepDecTop2); };
121  void addEventHypo(const HypoClassKey& key, const HypoCombPair hyp) { evtHyp_[key].push_back(hyp); };
123  void setNumberOfConsideredJets(const HypoClassKey& key, const unsigned int nJets) { nJetsConsidered_[key]=nJets; };
125  void setGenMatchSumPt(const std::vector<double>& val) {genMatchSumPt_=val;};
127  void setGenMatchSumDR(const std::vector<double>& val) {genMatchSumDR_=val;};
129  void setMvaMethod(const std::string& name) { mvaMethod_=name; };
131  void setMvaDiscriminators(const std::vector<double>& val) { mvaDisc_=val; };
133  void setFitChi2(const std::vector<double>& val) { fitChi2_=val; };
135  void setHitFitChi2(const std::vector<double>& val) { hitFitChi2_=val; };
137  void setFitProb(const std::vector<double>& val) { fitProb_=val; };
139  void setHitFitProb(const std::vector<double>& val) { hitFitProb_=val; };
141  void setHitFitMT(const std::vector<double>& val) { hitFitMT_=val; };
143  void setHitFitSigMT(const std::vector<double>& val) { hitFitSigMT_=val; };
144 
145  protected:
146 
148  std::pair<WDecay::LeptonType, WDecay::LeptonType> lepDecays_;
153  std::map<HypoClassKey, std::vector<HypoCombPair> > evtHyp_;
155  std::map<HypoClassKey, int> nJetsConsidered_;
156 
158  std::vector<double> fitChi2_;
159  std::vector<double> hitFitChi2_;
161  std::vector<double> fitProb_;
162  std::vector<double> hitFitProb_;
164  std::vector<double> hitFitMT_;
165  std::vector<double> hitFitSigMT_;
167  std::vector<double> genMatchSumPt_;
169  std::vector<double> genMatchSumDR_;
171  std::string mvaMethod_;
173  std::vector<double> mvaDisc_;
174 };
175 
176 #endif
void setMvaDiscriminators(const std::vector< double > &val)
set mva discriminant values of kMVADisc hypothesis
Definition: TtEvent.h:131
std::map< HypoClassKey, int > nJetsConsidered_
number of jets considered when building the hypotheses
Definition: TtEvent.h:155
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:62
HypoClassKey
supported classes of event hypotheses
Definition: TtEvent.h:31
void setHitFitSigMT(const std::vector< double > &val)
set fitted top mass uncertainty of kHitFit hypothesis
Definition: TtEvent.h:143
bool isHypoAvailable(const std::string &key, const unsigned &cmb=0) const
Definition: TtEvent.h:60
std::vector< double > genMatchSumDR_
result of gen match
Definition: TtEvent.h:169
TtEvent()
empty constructor
Definition: TtEvent.h:43
const edm::RefProd< TtGenEvent > & genEvent() const
get TtGenEvent
Definition: TtEvent.h:53
const math::XYZTLorentzVector * topPair() const
get combined 4-vector of top and topBar from the TtGenEvent
Definition: TtEvent.h:111
unsigned int numberOfAvailableHypos(const std::string &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:70
role_collection const & roles() const
get the roles
void setMvaMethod(const std::string &name)
set label of mva method for kMVADisc hypothesis
Definition: TtEvent.h:129
std::vector< double > hitFitChi2_
Definition: TtEvent.h:159
void setGenMatchSumDR(const std::vector< double > &val)
set sum dr of kGenMatch hypothesis
Definition: TtEvent.h:127
void setFitChi2(const std::vector< double > &val)
set chi2 of kKinFit hypothesis
Definition: TtEvent.h:133
bool isHypoClassAvailable(const std::string &key) const
check if hypothesis class &#39;key&#39; was added to the event structure
Definition: TtEvent.h:56
std::pair< WDecay::LeptonType, WDecay::LeptonType > lepDecays() const
get leptonic decay channels
Definition: TtEvent.h:48
unsigned int numberOfAvailableHypoClasses() const
return number of available hypothesis classes
Definition: TtEvent.h:68
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:80
void setLepDecays(const WDecay::LeptonType &lepDecTop1, const WDecay::LeptonType &lepDecTop2)
set leptonic decay channels
Definition: TtEvent.h:117
std::pair< reco::CompositeCandidate, std::vector< int > > HypoCombPair
pair of hypothesis and lepton jet combinatorics for a given hypothesis
Definition: TtEvent.h:33
std::vector< double > mvaDisc_
MVA discriminants.
Definition: TtEvent.h:173
std::pair< WDecay::LeptonType, WDecay::LeptonType > lepDecays_
leptonic decay channels
Definition: TtEvent.h:143
void printParticle(edm::LogInfo &log, const char *name, const reco::Candidate *cand) const
print pt, eta, phi, mass of a given candidate into an existing LogInfo
Definition: TtEvent.cc:51
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:66
U second(std::pair< T, U > const &p)
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:90
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
virtual ~TtEvent()
default destructor
Definition: TtEvent.h:45
double genMatchSumDR(const unsigned &cmb=0) const
return the sum dr of the generator match if available; -1 else
Definition: TtEvent.h:84
int numberOfConsideredJets(const std::string &key) const
return number of jets that were considered when building a given hypothesis
Definition: TtEvent.h:74
std::vector< double > fitProb_
result of kinematic fit
Definition: TtEvent.h:161
std::string mvaMethod_
label of the MVA method
Definition: TtEvent.h:171
void setGenMatchSumPt(const std::vector< double > &val)
set sum pt of kGenMatch hypothesis
Definition: TtEvent.h:125
double hitFitMT(const unsigned &cmb=0) const
return the hitfit top mass of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:98
void setHitFitMT(const std::vector< double > &val)
set fitted top mass of kHitFit hypothesis
Definition: TtEvent.h:141
void setHitFitChi2(const std::vector< double > &val)
set chi2 of kHitFit hypothesis
Definition: TtEvent.h:135
std::vector< double > genMatchSumPt_
result of gen match
Definition: TtEvent.h:167
HypoClassKey hypoClassKeyFromString(const std::string &label) const
return the corresponding enum value from a string
Definition: TtEvent.cc:18
std::vector< double > hitFitSigMT_
Definition: TtEvent.h:165
void setFitProb(const std::vector< double > &val)
set fit probability of kKinFit hypothesis
Definition: TtEvent.h:137
bool first
Definition: L1TdeRCT.cc:94
void setGenEvent(const edm::Handle< TtGenEvent > &evt)
set TtGenEvent
Definition: TtEvent.h:119
std::vector< double > hitFitMT_
result of hitfit
Definition: TtEvent.h:164
a lightweight map for selection type string label and enum value
Definition: TtEvent.h:37
double mvaDisc(const unsigned &cmb=0) const
return the mva discriminant value of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:88
std::vector< double > hitFitProb_
Definition: TtEvent.h:162
const reco::CompositeCandidate & eventHypo(const HypoClassKey &key, const unsigned &cmb=0) const
Definition: TtEvent.h:51
Base class to hold information for ttbar event interpretation.
Definition: TtEvent.h:27
bool isHypoClassAvailable(const HypoClassKey &key) const
check if hypothesis class &#39;key&#39; was added to the event structure
Definition: TtEvent.h:58
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:102
edm::RefProd< TtGenEvent > genEvt_
reference to TtGenEvent (has to be kept in the event!)
Definition: TtEvent.h:150
double hitFitChi2(const unsigned &cmb=0) const
return the hitfit chi2 of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:92
void setNumberOfConsideredJets(const HypoClassKey &key, const unsigned int nJets)
set number of jets considered when building a given hypothesis
Definition: TtEvent.h:123
double fitProb(const unsigned &cmb=0) const
return the fit probability of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:94
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:64
double genMatchSumPt(const unsigned &cmb=0) const
return the sum pt of the generator match if available; -1 else
Definition: TtEvent.h:82
std::vector< double > fitChi2_
result of kinematic fit
Definition: TtEvent.h:158
list key
Definition: combine.py:13
double hitFitProb(const unsigned &cmb=0) const
return the hitfit probability of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:96
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:107
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:86
int numberOfConsideredJets(const HypoClassKey &key) const
return number of jets that were considered when building a given hypothesis
Definition: TtEvent.h:76
std::map< HypoClassKey, std::vector< HypoCombPair > > evtHyp_
Definition: TtEvent.h:153
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:78
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:109
void setHitFitProb(const std::vector< double > &val)
set fit probability of kHitFit hypothesis
Definition: TtEvent.h:139
unsigned int numberOfAvailableHypos(const HypoClassKey &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:72
void addEventHypo(const HypoClassKey &key, const HypoCombPair hyp)
add new hypotheses
Definition: TtEvent.h:121
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:100