Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 //
0002 //
0003 
0004 #ifndef TopObjects_TtDilepEvtSolution_h
0005 #define TopObjects_TtDilepEvtSolution_h
0006 
0007 #include <vector>
0008 #include <string>
0009 
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "DataFormats/Common/interface/RefProd.h"
0012 #include "DataFormats/Common/interface/Ref.h"
0013 
0014 #include "DataFormats/Candidate/interface/Particle.h"
0015 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0016 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
0017 
0018 #include "DataFormats/PatCandidates/interface/Particle.h"
0019 #include "DataFormats/PatCandidates/interface/Electron.h"
0020 #include "DataFormats/PatCandidates/interface/Muon.h"
0021 #include "DataFormats/PatCandidates/interface/Tau.h"
0022 #include "DataFormats/PatCandidates/interface/Jet.h"
0023 #include "DataFormats/PatCandidates/interface/MET.h"
0024 
0025 class TtDilepEvtSolution {
0026   friend class TtFullLepKinSolver;
0027   friend class TtDilepEvtSolutionMaker;
0028   friend class TtDilepLRSignalSelObservables;
0029   friend class TtLRSignalSelCalc;
0030 
0031 public:
0032   TtDilepEvtSolution();
0033   virtual ~TtDilepEvtSolution();
0034 
0035   //-------------------------------------------
0036   // get calibrated base objects
0037   //-------------------------------------------
0038   pat::Jet getJetB() const;
0039   pat::Jet getJetBbar() const;
0040   pat::Electron getElectronp() const { return *elecp_; };
0041   pat::Electron getElectronm() const { return *elecm_; };
0042   pat::Muon getMuonp() const { return *muonp_; };
0043   pat::Muon getMuonm() const { return *muonm_; };
0044   pat::Tau getTaup() const { return *taup_; };
0045   pat::Tau getTaum() const { return *taum_; };
0046   pat::MET getMET() const { return *met_; };
0047 
0048   //-------------------------------------------
0049   // get the matched gen particles
0050   //-------------------------------------------
0051   const edm::RefProd<TtGenEvent>& getGenEvent() const { return theGenEvt_; };
0052   const reco::GenParticle* getGenT() const {
0053     if (!theGenEvt_)
0054       return nullptr;
0055     else
0056       return theGenEvt_->top();
0057   };
0058   const reco::GenParticle* getGenWp() const {
0059     if (!theGenEvt_)
0060       return nullptr;
0061     else
0062       return theGenEvt_->wPlus();
0063   };
0064   const reco::GenParticle* getGenB() const {
0065     if (!theGenEvt_)
0066       return nullptr;
0067     else
0068       return theGenEvt_->b();
0069   };
0070   const reco::GenParticle* getGenLepp() const {
0071     if (!theGenEvt_)
0072       return nullptr;
0073     else
0074       return theGenEvt_->leptonBar();
0075   };
0076   const reco::GenParticle* getGenN() const {
0077     if (!theGenEvt_)
0078       return nullptr;
0079     else
0080       return theGenEvt_->neutrino();
0081   };
0082   const reco::GenParticle* getGenTbar() const {
0083     if (!theGenEvt_)
0084       return nullptr;
0085     else
0086       return theGenEvt_->topBar();
0087   };
0088   const reco::GenParticle* getGenWm() const {
0089     if (!theGenEvt_)
0090       return nullptr;
0091     else
0092       return theGenEvt_->wMinus();
0093   };
0094   const reco::GenParticle* getGenBbar() const {
0095     if (!theGenEvt_)
0096       return nullptr;
0097     else
0098       return theGenEvt_->bBar();
0099   };
0100   const reco::GenParticle* getGenLepm() const {
0101     if (!theGenEvt_)
0102       return nullptr;
0103     else
0104       return theGenEvt_->lepton();
0105   };
0106   const reco::GenParticle* getGenNbar() const {
0107     if (!theGenEvt_)
0108       return nullptr;
0109     else
0110       return theGenEvt_->neutrinoBar();
0111   };
0112 
0113   //-------------------------------------------
0114   // get (un-)/calibrated reco objects
0115   //-------------------------------------------
0116   pat::Jet getRecJetB() const { return this->getJetB().correctedJet("RAW"); };
0117   pat::Jet getCalJetB() const { return this->getJetB(); };
0118   pat::Jet getRecJetBbar() const { return this->getJetBbar().correctedJet("RAW"); };
0119   pat::Jet getCalJetBbar() const { return this->getJetBbar(); };
0120 
0121   //-------------------------------------------
0122   // get info on the W decays
0123   //-------------------------------------------
0124   std::string getWpDecay() const { return wpDecay_; }
0125   std::string getWmDecay() const { return wmDecay_; }
0126 
0127   //-------------------------------------------
0128   // miscellaneous
0129   //-------------------------------------------
0130   double getJetResidual() const;
0131   double getLeptonResidual() const;
0132   double getFullResidual() const { return getJetResidual() + getLeptonResidual(); }
0133   bool getBestSol() const { return bestSol_; }
0134   double getRecTopMass() const { return topmass_; }
0135   double getRecWeightMax() const { return weightmax_; }
0136 
0137   //-------------------------------------------
0138   // returns the 4-vector of the positive
0139   // lepton, with the charge and the pdgId
0140   //-------------------------------------------
0141   reco::Particle getLeptPos() const;
0142 
0143   //-------------------------------------------
0144   // returns the 4-vector of the negative
0145   // lepton, with the charge and the pdgId
0146   //-------------------------------------------
0147   reco::Particle getLeptNeg() const;
0148 
0149   //-------------------------------------------
0150   // get info on the outcome of the signal
0151   // selection LR
0152   //-------------------------------------------
0153   double getLRSignalEvtObsVal(unsigned int) const;
0154   double getLRSignalEvtLRval() const { return lrSignalEvtLRval_; }
0155   double getLRSignalEvtProb() const { return lrSignalEvtProb_; }
0156 
0157 protected:
0158   //-------------------------------------------
0159   // set the generated event
0160   //-------------------------------------------
0161   void setGenEvt(const edm::Handle<TtGenEvent>&);
0162 
0163   //-------------------------------------------
0164   // set the basic objects
0165   //-------------------------------------------
0166   void setJetCorrectionScheme(int jetCorrScheme) { jetCorrScheme_ = jetCorrScheme; };
0167   void setB(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0168     jetB_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0169   };
0170   void setBbar(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0171     jetBbar_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0172   };
0173   void setMuonp(const edm::Handle<std::vector<pat::Muon> >& muon, int i) {
0174     muonp_ = edm::Ref<std::vector<pat::Muon> >(muon, i);
0175     wpDecay_ = "muon";
0176   };
0177   void setMuonm(const edm::Handle<std::vector<pat::Muon> >& muon, int i) {
0178     muonm_ = edm::Ref<std::vector<pat::Muon> >(muon, i);
0179     wmDecay_ = "muon";
0180   }
0181   void setTaup(const edm::Handle<std::vector<pat::Tau> >& tau, int i) {
0182     taup_ = edm::Ref<std::vector<pat::Tau> >(tau, i);
0183     wpDecay_ = "tau";
0184   }
0185   void setTaum(const edm::Handle<std::vector<pat::Tau> >& tau, int i) {
0186     taum_ = edm::Ref<std::vector<pat::Tau> >(tau, i);
0187     wmDecay_ = "tau";
0188   }
0189   void setElectronp(const edm::Handle<std::vector<pat::Electron> >& elec, int i) {
0190     elecp_ = edm::Ref<std::vector<pat::Electron> >(elec, i);
0191     wpDecay_ = "electron";
0192   };
0193   void setElectronm(const edm::Handle<std::vector<pat::Electron> >& elec, int i) {
0194     elecm_ = edm::Ref<std::vector<pat::Electron> >(elec, i);
0195     wmDecay_ = "electron";
0196   };
0197   void setMET(const edm::Handle<std::vector<pat::MET> >& met, int i) {
0198     met_ = edm::Ref<std::vector<pat::MET> >(met, i);
0199   };
0200 
0201   //-------------------------------------------
0202   // miscellaneous
0203   //-------------------------------------------
0204   void setBestSol(bool bs) { bestSol_ = bs; };
0205   void setRecTopMass(double mass) { topmass_ = mass; };
0206   void setRecWeightMax(double wgt) { weightmax_ = wgt; };
0207 
0208   //-------------------------------------------
0209   // set the outcome of the signal selection LR
0210   //-------------------------------------------
0211   void setLRSignalEvtObservables(const std::vector<std::pair<unsigned int, double> >&);
0212   void setLRSignalEvtLRval(double clr) { lrSignalEvtLRval_ = clr; };
0213   void setLRSignalEvtProb(double plr) { lrSignalEvtProb_ = plr; };
0214 
0215 private:
0216   //-------------------------------------------
0217   // particle content
0218   //-------------------------------------------
0219   edm::RefProd<TtGenEvent> theGenEvt_;
0220   edm::Ref<std::vector<pat::Electron> > elecp_, elecm_;
0221   edm::Ref<std::vector<pat::Muon> > muonp_, muonm_;
0222   edm::Ref<std::vector<pat::Tau> > taup_, taum_;
0223   edm::Ref<std::vector<pat::Jet> > jetB_, jetBbar_;
0224   edm::Ref<std::vector<pat::MET> > met_;
0225 
0226   //-------------------------------------------
0227   // miscellaneous
0228   //-------------------------------------------
0229   int jetCorrScheme_;
0230   std::string wpDecay_;
0231   std::string wmDecay_;
0232   bool bestSol_;
0233   double topmass_;
0234   double weightmax_;
0235 
0236   double lrSignalEvtLRval_, lrSignalEvtProb_;
0237   std::vector<std::pair<unsigned int, double> > lrSignalEvtVarVal_;
0238 };
0239 
0240 #endif