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  std::vector<int> jetLeptonCombination(const std::string& key, const unsigned& cmb=0) const { return jetLeptonCombination(hypoClassKeyFromString(key), cmb); };
76  std::vector<int> jetLeptonCombination(const HypoClassKey& key, const unsigned& cmb=0) const { return (evtHyp_.find(key)->second)[cmb].second; };
78  double genMatchSumPt(const unsigned& cmb=0) const { return (cmb<genMatchSumPt_.size() ? genMatchSumPt_[cmb] : -1.); };
80  double genMatchSumDR(const unsigned& cmb=0) const { return (cmb<genMatchSumDR_.size() ? genMatchSumDR_[cmb] : -1.); };
82  std::string mvaMethod() const { return mvaMethod_; }
84  double mvaDisc(const unsigned& cmb=0) const { return (cmb<mvaDisc_.size() ? mvaDisc_[cmb] : -1.); }
86  double fitChi2(const unsigned& cmb=0) const { return (cmb<fitChi2_.size() ? fitChi2_[cmb] : -1.); }
88  double fitProb(const unsigned& cmb=0) const { return (cmb<fitProb_.size() ? fitProb_[cmb] : -1.); }
90  int correspondingHypo(const std::string& key1, const unsigned& hyp1, const std::string& key2) const { return correspondingHypo(hypoClassKeyFromString(key1), hyp1, hypoClassKeyFromString(key2) ); };
92  int correspondingHypo(const HypoClassKey& key1, const unsigned& hyp1, const HypoClassKey& key2) const;
93 
95  void printParticle(edm::LogInfo &log, const char* name, const reco::Candidate* cand) const;
96 
98  void setLepDecays(const WDecay::LeptonType& lepDecTop1, const WDecay::LeptonType& lepDecTop2) { lepDecays_=std::make_pair(lepDecTop1, lepDecTop2); };
102  void addEventHypo(const HypoClassKey& key, HypoCombPair hyp) { evtHyp_[key].push_back(hyp); };
104  void setGenMatchSumPt(const std::vector<double>& val) {genMatchSumPt_=val;};
106  void setGenMatchSumDR(const std::vector<double>& val) {genMatchSumDR_=val;};
108  void setMvaMethod(const std::string& name) { mvaMethod_=name; };
110  void setMvaDiscriminators(const std::vector<double>& val) { mvaDisc_=val; };
112  void setFitChi2(const std::vector<double>& val) { fitChi2_=val; };
114  void setFitProb(const std::vector<double>& val) { fitProb_=val; };
115 
116  protected:
117 
119  std::pair<WDecay::LeptonType, WDecay::LeptonType> lepDecays_;
124  std::map<HypoClassKey, std::vector<HypoCombPair> > evtHyp_;
125 
127  std::vector<double> fitChi2_;
129  std::vector<double> fitProb_;
131  std::vector<double> genMatchSumPt_;
133  std::vector<double> genMatchSumDR_;
135  std::string mvaMethod_;
137  std::vector<double> mvaDisc_;
138 };
139 
140 #endif
void setMvaDiscriminators(const std::vector< double > &val)
set mva discriminant values of kMVADisc hypothesis
Definition: TtEvent.h:110
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
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:133
TtEvent()
empty constructor
Definition: TtEvent.h:43
const edm::RefProd< TtGenEvent > & genEvent() const
get TtGenEvent
Definition: TtEvent.h:53
const std::string & label
Definition: MVAComputer.cc:186
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:108
void setGenMatchSumDR(const std::vector< double > &val)
set sum dr of kGenMatch hypothesis
Definition: TtEvent.h:106
void setFitChi2(const std::vector< double > &val)
set chi2 of kKinFit hypothesis
Definition: TtEvent.h:112
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:76
void setLepDecays(const WDecay::LeptonType &lepDecTop1, const WDecay::LeptonType &lepDecTop2)
set leptonic decay channels
Definition: TtEvent.h:98
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:137
std::pair< WDecay::LeptonType, WDecay::LeptonType > lepDecays_
leptonic decay channels
Definition: TtEvent.h:114
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:50
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 kinemtaic fit of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:86
void addEventHypo(const HypoClassKey &key, HypoCombPair hyp)
add new hypotheses
Definition: TtEvent.h:102
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:80
std::vector< double > fitProb_
result of kinematic fit
Definition: TtEvent.h:129
std::string mvaMethod_
label of the MVA method
Definition: TtEvent.h:135
void setGenMatchSumPt(const std::vector< double > &val)
set sum pt of kGenMatch hypothesis
Definition: TtEvent.h:104
std::vector< double > genMatchSumPt_
result of gen match
Definition: TtEvent.h:131
HypoClassKey hypoClassKeyFromString(const std::string &label) const
return the corresponding enum value from a string
Definition: TtEvent.cc:18
void setFitProb(const std::vector< double > &val)
set fit probability of kKinFit hypothesis
Definition: TtEvent.h:114
bool first
Definition: L1TdeRCT.cc:79
void setGenEvent(const edm::Handle< TtGenEvent > &evt)
set TtGenEvent
Definition: TtEvent.h:100
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:84
const reco::CompositeCandidate & eventHypo(const HypoClassKey &key, const unsigned &cmb=0) const
Definition: TtEvent.h:51
Log< T >::type log(const T &t)
Definition: Log.h:22
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:90
edm::RefProd< TtGenEvent > genEvt_
reference to TtGenEvent (has to be kept in the event!)
Definition: TtEvent.h:121
double fitProb(const unsigned &cmb=0) const
return the fit probability of hypothesis &#39;cmb&#39; if available; -1 else
Definition: TtEvent.h:88
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:78
std::vector< double > fitChi2_
result of kinematic fit
Definition: TtEvent.h:127
list key
Definition: combine.py:13
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:82
std::map< HypoClassKey, std::vector< HypoCombPair > > evtHyp_
Definition: TtEvent.h:124
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:74
unsigned int numberOfAvailableHypos(const HypoClassKey &key) const
return number of available hypotheses within a given hypothesis class
Definition: TtEvent.h:72