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_StEvtSolution_h
0005 #define TopObjects_StEvtSolution_h
0006 
0007 #include <vector>
0008 
0009 #include "DataFormats/Common/interface/Handle.h"
0010 #include "DataFormats/Common/interface/RefProd.h"
0011 #include "DataFormats/Common/interface/Ref.h"
0012 
0013 #include "DataFormats/PatCandidates/interface/Particle.h"
0014 #include "DataFormats/PatCandidates/interface/Electron.h"
0015 #include "DataFormats/PatCandidates/interface/Muon.h"
0016 #include "DataFormats/PatCandidates/interface/Jet.h"
0017 #include "DataFormats/PatCandidates/interface/MET.h"
0018 
0019 #include "AnalysisDataFormats/TopObjects/interface/StGenEvent.h"
0020 
0021 class StEvtSolution {
0022   friend class StEvtSolutionMaker;
0023   friend class StKinFitter;
0024 
0025 public:
0026   StEvtSolution();
0027   virtual ~StEvtSolution();
0028 
0029   //-------------------------------------------
0030   // get calibrated base objects
0031   //-------------------------------------------
0032   pat::Jet getBottom() const;
0033   pat::Jet getLight() const;
0034   pat::Muon getMuon() const { return *muon_; };
0035   pat::Electron getElectron() const { return *electron_; };
0036   pat::MET getNeutrino() const { return *neutrino_; };
0037   reco::Particle getLepW() const;
0038   reco::Particle getLept() const;
0039 
0040   //-------------------------------------------
0041   // get the matched gen particles
0042   //-------------------------------------------
0043   const edm::RefProd<StGenEvent>& getGenEvent() const { return theGenEvt_; };
0044   const reco::GenParticle* getGenBottom() const;
0045   //const reco::GenParticle * getGenLight() const; // not implemented yet
0046   const reco::GenParticle* getGenLepton() const;
0047   const reco::GenParticle* getGenNeutrino() const;
0048   const reco::GenParticle* getGenLepW() const;
0049   const reco::GenParticle* getGenLept() const;
0050 
0051   //-------------------------------------------
0052   // get uncalibrated reco objects
0053   //-------------------------------------------
0054   pat::Jet getRecBottom() const { return this->getBottom().correctedJet("RAW"); };
0055   pat::Jet getRecLight() const { return this->getLight().correctedJet("RAW"); };
0056   pat::Muon getRecMuon() const { return this->getMuon(); };              // redundant
0057   pat::Electron getRecElectron() const { return this->getElectron(); };  // redundant
0058   pat::MET getRecNeutrino() const { return this->getNeutrino(); };       // redundant
0059   reco::Particle getRecLepW() const { return this->getLepW(); };         // redundant
0060   reco::Particle getRecLept() const;
0061 
0062   //-------------------------------------------
0063   // get objects from kinematic fit
0064   //-------------------------------------------
0065   pat::Particle getFitBottom() const { return (!fitBottom_.empty() ? fitBottom_.front() : pat::Particle()); };
0066   pat::Particle getFitLight() const { return (!fitLight_.empty() ? fitLight_.front() : pat::Particle()); };
0067   pat::Particle getFitLepton() const { return (!fitLepton_.empty() ? fitLepton_.front() : pat::Particle()); };
0068   pat::Particle getFitNeutrino() const { return (!fitNeutrino_.empty() ? fitNeutrino_.front() : pat::Particle()); };
0069   reco::Particle getFitLepW() const;
0070   reco::Particle getFitLept() const;
0071 
0072   //-------------------------------------------
0073   // get info on the selected decay
0074   //-------------------------------------------
0075   std::string getDecay() const { return decay_; }
0076 
0077   //-------------------------------------------
0078   // get other event info
0079   //-------------------------------------------
0080   std::vector<double> getScanValues() const { return scanValues_; }
0081   double getChi2Prob() const { return chi2Prob_; }
0082   double getPtrueCombExist() const { return pTrueCombExist_; }
0083   double getPtrueBJetSel() const { return pTrueBJetSel_; }
0084   double getPtrueBhadrSel() const { return pTrueBhadrSel_; }
0085   double getPtrueJetComb() const { return pTrueJetComb_; }
0086   double getSignalPur() const { return signalPur_; }
0087   double getSignalLRTot() const { return signalLRTot_; }
0088   double getSumDeltaRjp() const { return sumDeltaRjp_; }
0089   double getDeltaRB() const { return deltaRB_; }
0090   double getDeltaRL() const { return deltaRL_; }
0091   int getChangeBL() const { return changeBL_; }
0092   bool getBestSol() const { return bestSol_; }
0093 
0094 protected:
0095   //-------------------------------------------
0096   // set the generated event
0097   //-------------------------------------------
0098   void setGenEvt(const edm::Handle<StGenEvent>&);
0099 
0100   //-------------------------------------------
0101   // set the basic objects
0102   //-------------------------------------------
0103   void setJetCorrectionScheme(int scheme) { jetCorrScheme_ = scheme; };
0104   void setBottom(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0105     bottom_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0106   };
0107   void setLight(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0108     light_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0109   };
0110   void setMuon(const edm::Handle<std::vector<pat::Muon> >& muon, int i) {
0111     muon_ = edm::Ref<std::vector<pat::Muon> >(muon, i);
0112     decay_ = "muon";
0113   };
0114   void setElectron(const edm::Handle<std::vector<pat::Electron> >& elec, int i) {
0115     electron_ = edm::Ref<std::vector<pat::Electron> >(elec, i);
0116     decay_ = "electron";
0117   };
0118   void setNeutrino(const edm::Handle<std::vector<pat::MET> >& met, int i) {
0119     neutrino_ = edm::Ref<std::vector<pat::MET> >(met, i);
0120   };
0121 
0122   //-------------------------------------------
0123   // set the fitted objects
0124   //-------------------------------------------
0125   void setFitBottom(const pat::Particle& part) {
0126     fitBottom_.clear();
0127     fitBottom_.push_back(part);
0128   };
0129   void setFitLight(const pat::Particle& part) {
0130     fitLight_.clear();
0131     fitLight_.push_back(part);
0132   };
0133   void setFitLepton(const pat::Particle& part) {
0134     fitLepton_.clear();
0135     fitLepton_.push_back(part);
0136   };
0137   void setFitNeutrino(const pat::Particle& part) {
0138     fitNeutrino_.clear();
0139     fitNeutrino_.push_back(part);
0140   };
0141 
0142   //-------------------------------------------
0143   // set other info on the event
0144   //-------------------------------------------
0145   void setChi2Prob(double prob) { chi2Prob_ = prob; };
0146   void setScanValues(const std::vector<double>&);
0147   void setPtrueCombExist(double pce) { pTrueCombExist_ = pce; };
0148   void setPtrueBJetSel(double pbs) { pTrueBJetSel_ = pbs; };
0149   void setPtrueBhadrSel(double pbh) { pTrueBhadrSel_ = pbh; };
0150   void setPtrueJetComb(double pt) { pTrueJetComb_ = pt; };
0151   void setSignalPurity(double pur) { signalPur_ = pur; };
0152   void setSignalLRTot(double lrt) { signalLRTot_ = lrt; };
0153   void setSumDeltaRjp(double sdr) { sumDeltaRjp_ = sdr; };
0154   void setDeltaRB(double adr) { deltaRB_ = adr; };
0155   void setDeltaRL(double adr) { deltaRL_ = adr; };
0156   void setChangeBL(int bl) { changeBL_ = bl; };
0157   void setBestSol(bool bs) { bestSol_ = bs; };
0158 
0159 private:
0160   //-------------------------------------------
0161   // particle content
0162   //-------------------------------------------
0163   edm::RefProd<StGenEvent> theGenEvt_;
0164   edm::Ref<std::vector<pat::Jet> > bottom_, light_;
0165   edm::Ref<std::vector<pat::Muon> > muon_;
0166   edm::Ref<std::vector<pat::Electron> > electron_;
0167   edm::Ref<std::vector<pat::MET> > neutrino_;
0168   std::vector<pat::Particle> fitBottom_, fitLight_, fitLepton_, fitNeutrino_;
0169 
0170   //-------------------------------------------
0171   // miscellaneous
0172   //-------------------------------------------
0173   std::string decay_;
0174   int jetCorrScheme_;
0175   double chi2Prob_;
0176   std::vector<double> scanValues_;
0177   double pTrueCombExist_, pTrueBJetSel_, pTrueBhadrSel_, pTrueJetComb_;
0178   double signalPur_, signalLRTot_;
0179   double sumDeltaRjp_, deltaRB_, deltaRL_;
0180   int changeBL_;
0181   bool bestSol_;
0182   //double jetMatchPur_;
0183 };
0184 
0185 #endif