Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:31

0001 #ifndef TopObjects_TtEvent_h
0002 #define TopObjects_TtEvent_h
0003 
0004 #include <vector>
0005 #include <string>
0006 
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "DataFormats/Common/interface/RefProd.h"
0009 #include "DataFormats/Candidate/interface/CompositeCandidate.h"
0010 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
0011 
0012 /**
0013    \class   TtEvent TtEvent.h "AnalysisDataFormats/TopObjects/interface/TtEvent.h"
0014 
0015    \brief   Base class to hold information for ttbar event interpretation
0016 
0017    The structure holds information for ttbar event interpretation.
0018    All event hypotheses of different classes (user defined during
0019    production) and a reference to the TtGenEvent (if available). It 
0020    provides access and administration.
0021 */
0022 
0023 class TtEvent {
0024 public:
0025   /// supported classes of event hypotheses
0026   enum HypoClassKey {
0027     kGeom,
0028     kWMassMaxSumPt,
0029     kMaxSumPtWMass,
0030     kGenMatch,
0031     kMVADisc,
0032     kKinFit,
0033     kKinSolution,
0034     kWMassDeltaTopMass,
0035     kHitFit
0036   };
0037   /// pair of hypothesis and lepton jet combinatorics for a given hypothesis
0038   typedef std::pair<reco::CompositeCandidate, std::vector<int> > HypoCombPair;
0039   /// a lightweight map for selection type string label and enum value
0040   struct HypoClassKeyStringToEnum {
0041     const char* label;
0042     HypoClassKey value;
0043   };
0044 
0045 protected:
0046   /// return the corresponding enum value from a string
0047   HypoClassKey hypoClassKeyFromString(const std::string& label) const;
0048 
0049 public:
0050   /// empty constructor
0051   TtEvent(){};
0052   /// default destructor
0053   virtual ~TtEvent(){};
0054 
0055   /// get leptonic decay channels
0056   std::pair<WDecay::LeptonType, WDecay::LeptonType> lepDecays() const { return lepDecays_; }
0057   /// get event hypothesis; there can be more hypotheses of a certain
0058   /// class (sorted by quality); per default the best hypothesis is returned
0059   const reco::CompositeCandidate& eventHypo(const HypoClassKey& key, const unsigned& cmb = 0) const {
0060     return (evtHyp_.find(key)->second)[cmb].first;
0061   };
0062   /// get TtGenEvent
0063   const edm::RefProd<TtGenEvent>& genEvent() const { return genEvt_; };
0064 
0065   /// check if hypothesis class 'key' was added to the event structure
0066   bool isHypoClassAvailable(const std::string& key) const { return isHypoClassAvailable(hypoClassKeyFromString(key)); };
0067   /// check if hypothesis class 'key' was added to the event structure
0068   bool isHypoClassAvailable(const HypoClassKey& key) const { return (evtHyp_.find(key) != evtHyp_.end()); };
0069   // check if hypothesis 'cmb' is available within the hypothesis class
0070   bool isHypoAvailable(const std::string& key, const unsigned& cmb = 0) const {
0071     return isHypoAvailable(hypoClassKeyFromString(key), cmb);
0072   };
0073   /// check if hypothesis 'cmb' is available within the hypothesis class
0074   bool isHypoAvailable(const HypoClassKey& key, const unsigned& cmb = 0) const {
0075     return isHypoClassAvailable(key) ? (cmb < evtHyp_.find(key)->second.size()) : false;
0076   };
0077   /// check if hypothesis 'cmb' within the hypothesis class was valid; if not it lead to an empty CompositeCandidate
0078   bool isHypoValid(const std::string& key, const unsigned& cmb = 0) const {
0079     return isHypoValid(hypoClassKeyFromString(key), cmb);
0080   };
0081   /// check if hypothesis 'cmb' within the hypothesis class was valid; if not it lead to an empty CompositeCandidate
0082   bool isHypoValid(const HypoClassKey& key, const unsigned& cmb = 0) const {
0083     return isHypoAvailable(key, cmb) ? !eventHypo(key, cmb).roles().empty() : false;
0084   };
0085   /// return number of available hypothesis classes
0086   unsigned int numberOfAvailableHypoClasses() const { return evtHyp_.size(); };
0087   /// return number of available hypotheses within a given hypothesis class
0088   unsigned int numberOfAvailableHypos(const std::string& key) const {
0089     return numberOfAvailableHypos(hypoClassKeyFromString(key));
0090   };
0091   /// return number of available hypotheses within a given hypothesis class
0092   unsigned int numberOfAvailableHypos(const HypoClassKey& key) const {
0093     return isHypoAvailable(key) ? evtHyp_.find(key)->second.size() : 0;
0094   };
0095   /// return number of jets that were considered when building a given hypothesis
0096   int numberOfConsideredJets(const std::string& key) const {
0097     return numberOfConsideredJets(hypoClassKeyFromString(key));
0098   };
0099   /// return number of jets that were considered when building a given hypothesis
0100   int numberOfConsideredJets(const HypoClassKey& key) const {
0101     return (isHypoAvailable(key) ? nJetsConsidered_.find(key)->second : -1);
0102   };
0103   /// return the vector of jet lepton combinatorics for a given hypothesis and class
0104   std::vector<int> jetLeptonCombination(const std::string& key, const unsigned& cmb = 0) const {
0105     return jetLeptonCombination(hypoClassKeyFromString(key), cmb);
0106   };
0107   /// return the vector of jet lepton combinatorics for a given hypothesis and class
0108   std::vector<int> jetLeptonCombination(const HypoClassKey& key, const unsigned& cmb = 0) const {
0109     return (evtHyp_.find(key)->second)[cmb].second;
0110   };
0111   /// return the sum pt of the generator match if available; -1 else
0112   double genMatchSumPt(const unsigned& cmb = 0) const {
0113     return (cmb < genMatchSumPt_.size() ? genMatchSumPt_[cmb] : -1.);
0114   };
0115   /// return the sum dr of the generator match if available; -1 else
0116   double genMatchSumDR(const unsigned& cmb = 0) const {
0117     return (cmb < genMatchSumDR_.size() ? genMatchSumDR_[cmb] : -1.);
0118   };
0119   /// return the label of the mva method in use for the jet parton association (if kMVADisc is not available the string is empty)
0120   std::string mvaMethod() const { return mvaMethod_; }
0121   /// return the mva discriminant value of hypothesis 'cmb' if available; -1 else
0122   double mvaDisc(const unsigned& cmb = 0) const { return (cmb < mvaDisc_.size() ? mvaDisc_[cmb] : -1.); }
0123   /// return the chi2 of the kinematic fit of hypothesis 'cmb' if available; -1 else
0124   double fitChi2(const unsigned& cmb = 0) const { return (cmb < fitChi2_.size() ? fitChi2_[cmb] : -1.); }
0125   /// return the hitfit chi2 of hypothesis 'cmb' if available; -1 else
0126   double hitFitChi2(const unsigned& cmb = 0) const { return (cmb < hitFitChi2_.size() ? hitFitChi2_[cmb] : -1.); }
0127   /// return the fit probability of hypothesis 'cmb' if available; -1 else
0128   double fitProb(const unsigned& cmb = 0) const { return (cmb < fitProb_.size() ? fitProb_[cmb] : -1.); }
0129   /// return the hitfit probability of hypothesis 'cmb' if available; -1 else
0130   double hitFitProb(const unsigned& cmb = 0) const { return (cmb < hitFitProb_.size() ? hitFitProb_[cmb] : -1.); }
0131   /// return the hitfit top mass of hypothesis 'cmb' if available; -1 else
0132   double hitFitMT(const unsigned& cmb = 0) const { return (cmb < hitFitMT_.size() ? hitFitMT_[cmb] : -1.); }
0133   /// return the hitfit top mass uncertainty of hypothesis 'cmb' if available; -1 else
0134   double hitFitSigMT(const unsigned& cmb = 0) const { return (cmb < hitFitSigMT_.size() ? hitFitSigMT_[cmb] : -1.); }
0135   /// return the hypothesis in hypothesis class 'key2', which corresponds to hypothesis 'hyp1' in hypothesis class 'key1'
0136   int correspondingHypo(const std::string& key1, const unsigned& hyp1, const std::string& key2) const {
0137     return correspondingHypo(hypoClassKeyFromString(key1), hyp1, hypoClassKeyFromString(key2));
0138   };
0139   /// return the hypothesis in hypothesis class 'key2', which corresponds to hypothesis 'hyp1' in hypothesis class 'key1'
0140   int correspondingHypo(const HypoClassKey& key1, const unsigned& hyp1, const HypoClassKey& key2) const;
0141 
0142   /// get combined 4-vector of top and topBar of the given hypothesis
0143   const reco::Candidate* topPair(const std::string& key, const unsigned& cmb = 0) const {
0144     return topPair(hypoClassKeyFromString(key), cmb);
0145   };
0146   /// get combined 4-vector of top and topBar of the given hypothesis
0147   const reco::Candidate* topPair(const HypoClassKey& key, const unsigned& cmb = 0) const {
0148     return !isHypoValid(key, cmb) ? nullptr : (reco::Candidate*)&eventHypo(key, cmb);
0149   };
0150   /// get combined 4-vector of top and topBar from the TtGenEvent
0151   const math::XYZTLorentzVector* topPair() const { return (!genEvt_ ? nullptr : this->genEvent()->topPair()); };
0152 
0153   /// set leptonic decay channels
0154   void setLepDecays(const WDecay::LeptonType& lepDecTop1, const WDecay::LeptonType& lepDecTop2) {
0155     lepDecays_ = std::make_pair(lepDecTop1, lepDecTop2);
0156   };
0157   /// set TtGenEvent
0158   void setGenEvent(const edm::Handle<TtGenEvent>& evt) { genEvt_ = edm::RefProd<TtGenEvent>(evt); };
0159   /// add new hypotheses
0160   void addEventHypo(const HypoClassKey& key, const HypoCombPair& hyp) { evtHyp_[key].push_back(hyp); };
0161   /// set number of jets considered when building a given hypothesis
0162   void setNumberOfConsideredJets(const HypoClassKey& key, const unsigned int nJets) { nJetsConsidered_[key] = nJets; };
0163   /// set sum pt of kGenMatch hypothesis
0164   void setGenMatchSumPt(const std::vector<double>& val) { genMatchSumPt_ = val; };
0165   /// set sum dr of kGenMatch hypothesis
0166   void setGenMatchSumDR(const std::vector<double>& val) { genMatchSumDR_ = val; };
0167   /// set label of mva method for kMVADisc hypothesis
0168   void setMvaMethod(const std::string& name) { mvaMethod_ = name; };
0169   /// set mva discriminant values of kMVADisc hypothesis
0170   void setMvaDiscriminators(const std::vector<double>& val) { mvaDisc_ = val; };
0171   /// set chi2 of kKinFit hypothesis
0172   void setFitChi2(const std::vector<double>& val) { fitChi2_ = val; };
0173   /// set chi2 of kHitFit hypothesis
0174   void setHitFitChi2(const std::vector<double>& val) { hitFitChi2_ = val; };
0175   /// set fit probability of kKinFit hypothesis
0176   void setFitProb(const std::vector<double>& val) { fitProb_ = val; };
0177   /// set fit probability of kHitFit hypothesis
0178   void setHitFitProb(const std::vector<double>& val) { hitFitProb_ = val; };
0179   /// set fitted top mass of kHitFit hypothesis
0180   void setHitFitMT(const std::vector<double>& val) { hitFitMT_ = val; };
0181   /// set fitted top mass uncertainty of kHitFit hypothesis
0182   void setHitFitSigMT(const std::vector<double>& val) { hitFitSigMT_ = val; };
0183 
0184 protected:
0185   /// leptonic decay channels
0186   std::pair<WDecay::LeptonType, WDecay::LeptonType> lepDecays_;
0187   /// reference to TtGenEvent (has to be kept in the event!)
0188   edm::RefProd<TtGenEvent> genEvt_;
0189   /// map of hypotheses; for each HypoClassKey a vector of
0190   /// hypothesis and their lepton jet combinatorics are kept
0191   std::map<HypoClassKey, std::vector<HypoCombPair> > evtHyp_;
0192   /// number of jets considered when building the hypotheses
0193   std::map<HypoClassKey, int> nJetsConsidered_;
0194 
0195   /// result of kinematic fit
0196   std::vector<double> fitChi2_;
0197   std::vector<double> hitFitChi2_;
0198   /// result of kinematic fit
0199   std::vector<double> fitProb_;
0200   std::vector<double> hitFitProb_;
0201   /// result of hitfit
0202   std::vector<double> hitFitMT_;
0203   std::vector<double> hitFitSigMT_;
0204   /// result of gen match
0205   std::vector<double> genMatchSumPt_;
0206   /// result of gen match
0207   std::vector<double> genMatchSumDR_;
0208   /// label of the MVA method
0209   std::string mvaMethod_;
0210   /// MVA discriminants
0211   std::vector<double> mvaDisc_;
0212 };
0213 
0214 #endif