CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/AnalysisDataFormats/TopObjects/interface/StEvtSolution.h

Go to the documentation of this file.
00001 //
00002 // $Id: StEvtSolution.h,v 1.13 2008/12/18 21:20:10 rwolf Exp $
00003 //
00004 
00005 #ifndef TopObjects_StEvtSolution_h
00006 #define TopObjects_StEvtSolution_h
00007 
00008 #include <vector>
00009 
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "DataFormats/Common/interface/RefProd.h"
00012 #include "DataFormats/Common/interface/Ref.h"
00013 
00014 #include "DataFormats/PatCandidates/interface/Particle.h"
00015 #include "DataFormats/PatCandidates/interface/Electron.h"
00016 #include "DataFormats/PatCandidates/interface/Muon.h"
00017 #include "DataFormats/PatCandidates/interface/Jet.h"
00018 #include "DataFormats/PatCandidates/interface/MET.h"
00019 
00020 #include "AnalysisDataFormats/TopObjects/interface/StGenEvent.h"
00021 
00022 class StEvtSolution {
00023 
00024   friend class StEvtSolutionMaker;
00025   friend class StKinFitter;
00026   
00027  public:
00028 
00029   StEvtSolution();
00030   virtual ~StEvtSolution();
00031   
00032   //-------------------------------------------
00033   // get calibrated base objects 
00034   //-------------------------------------------
00035   pat::Jet       getBottom()   const;
00036   pat::Jet       getLight()    const;
00037   pat::Muon      getMuon()     const { return *muon_; };
00038   pat::Electron  getElectron() const { return *electron_; };
00039   pat::MET       getNeutrino() const { return *neutrino_; };
00040   reco::Particle getLepW()     const;  
00041   reco::Particle getLept()     const;
00042 
00043   //-------------------------------------------
00044   // get the matched gen particles
00045   //-------------------------------------------
00046   const edm::RefProd<StGenEvent>& getGenEvent() const { return theGenEvt_; };
00047   const reco::GenParticle * getGenBottom()   const;
00048   //const reco::GenParticle * getGenLight() const; // not implemented yet
00049   const reco::GenParticle * getGenLepton()   const;
00050   const reco::GenParticle * getGenNeutrino() const;
00051   const reco::GenParticle * getGenLepW()     const;
00052   const reco::GenParticle * getGenLept()     const;
00053 
00054   //-------------------------------------------
00055   // get uncalibrated reco objects
00056   //-------------------------------------------
00057   pat::Jet       getRecBottom()   const { return this->getBottom().correctedJet("RAW"); };
00058   pat::Jet       getRecLight()    const { return this->getLight ().correctedJet("RAW"); };
00059   pat::Muon      getRecMuon()     const { return this->getMuon(); };     // redundant
00060   pat::Electron  getRecElectron() const { return this->getElectron(); }; // redundant
00061   pat::MET       getRecNeutrino() const { return this->getNeutrino(); }; // redundant
00062   reco::Particle getRecLepW()     const { return this->getLepW(); };     // redundant
00063   reco::Particle getRecLept()     const;
00064 
00065   //-------------------------------------------
00066   // get objects from kinematic fit
00067   //-------------------------------------------
00068   pat::Particle getFitBottom() const { return (fitBottom_.size()>0 ? fitBottom_.front() : pat::Particle()); };
00069   pat::Particle getFitLight() const { return (fitLight_.size()>0 ? fitLight_.front() : pat::Particle()); };
00070   pat::Particle getFitLepton() const { return (fitLepton_.size()>0 ? fitLepton_.front() : pat::Particle()); };
00071   pat::Particle getFitNeutrino() const { return (fitNeutrino_.size()>0 ? fitNeutrino_.front() : pat::Particle()); };
00072   reco::Particle getFitLepW() const;
00073   reco::Particle getFitLept() const;
00074 
00075   //-------------------------------------------
00076   // get info on the selected decay
00077   //-------------------------------------------
00078   std::string getDecay() const { return decay_; }
00079 
00080   //-------------------------------------------
00081   // get other event info
00082   //-------------------------------------------
00083   std::vector<double> getScanValues() const { return scanValues_; }
00084   double getChi2Prob()       const { return chi2Prob_; }
00085   double getPtrueCombExist() const { return pTrueCombExist_; }
00086   double getPtrueBJetSel()  const { return pTrueBJetSel_; }
00087   double getPtrueBhadrSel() const { return pTrueBhadrSel_; }
00088   double getPtrueJetComb() const { return pTrueJetComb_; }
00089   double getSignalPur()   const { return signalPur_; }
00090   double getSignalLRTot() const { return signalLRTot_; }
00091   double getSumDeltaRjp() const { return sumDeltaRjp_; }
00092   double getDeltaRB() const { return deltaRB_; }
00093   double getDeltaRL() const { return deltaRL_; }
00094   int  getChangeBL() const { return changeBL_; }
00095   bool getBestSol() const { return bestSol_; }
00096  
00097  protected:         
00098 
00099   //-------------------------------------------  
00100   // set the generated event
00101   //-------------------------------------------
00102   void setGenEvt(const edm::Handle<StGenEvent> &);
00103 
00104   //-------------------------------------------
00105   // set the basic objects
00106   //-------------------------------------------
00107   void setJetCorrectionScheme(int scheme) { jetCorrScheme_ = scheme;};
00108   void setBottom(const edm::Handle<std::vector<pat::Jet > >& jet, int i) 
00109   { bottom_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00110   void setLight (const edm::Handle<std::vector<pat::Jet > >& jet, int i) 
00111   { light_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00112   void setMuon  (const edm::Handle<std::vector<pat::Muon> >& muon, int i) 
00113   { muon_ = edm::Ref<std::vector<pat::Muon> >(muon, i); decay_ = "muon"; };
00114   void setElectron(const edm::Handle<std::vector<pat::Electron> >& elec, int i) 
00115   { electron_ = edm::Ref<std::vector<pat::Electron> >(elec, i); decay_ = "electron"; };
00116   void setNeutrino(const edm::Handle<std::vector<pat::MET> >& met, int i)
00117   { neutrino_ = edm::Ref<std::vector<pat::MET> >(met, i); };
00118 
00119   //-------------------------------------------
00120   // set the fitted objects 
00121   //-------------------------------------------
00122   void setFitBottom(const pat::Particle& part) { fitBottom_.clear(); fitBottom_.push_back(part); };
00123   void setFitLight (const pat::Particle& part) { fitLight_.clear(); fitLight_.push_back(part); };
00124   void setFitLepton(const pat::Particle& part) { fitLepton_.clear(); fitLepton_.push_back(part); };
00125   void setFitNeutrino(const pat::Particle& part) { fitNeutrino_.clear(); fitNeutrino_.push_back(part); };
00126 
00127   //-------------------------------------------
00128   // set other info on the event
00129   //-------------------------------------------
00130   void setChi2Prob(double prob){ chi2Prob_ = prob; };
00131   void setScanValues(const std::vector<double>&);
00132   void setPtrueCombExist(double pce){ pTrueCombExist_ = pce; };
00133   void setPtrueBJetSel (double pbs) { pTrueBJetSel_   = pbs; };
00134   void setPtrueBhadrSel(double pbh) { pTrueBhadrSel_  = pbh; };
00135   void setPtrueJetComb (double pt)  { pTrueJetComb_   = pt;  };
00136   void setSignalPurity (double pur) { signalPur_ = pur; };
00137   void setSignalLRTot(double lrt){ signalLRTot_ = lrt; };
00138   void setSumDeltaRjp(double sdr){ sumDeltaRjp_ = sdr; };
00139   void setDeltaRB(double adr) { deltaRB_ = adr; };
00140   void setDeltaRL(double adr) { deltaRL_ = adr; };
00141   void setChangeBL(int bl) { changeBL_ = bl; };
00142   void setBestSol(bool bs) { bestSol_  = bs; };
00143   
00144  private:
00145 
00146   //-------------------------------------------  
00147   // particle content
00148   //-------------------------------------------
00149   edm::RefProd<StGenEvent> theGenEvt_;
00150   edm::Ref<std::vector<pat::Jet> >  bottom_, light_;
00151   edm::Ref<std::vector<pat::Muon> > muon_;
00152   edm::Ref<std::vector<pat::Electron> > electron_;
00153   edm::Ref<std::vector<pat::MET> > neutrino_;
00154   std::vector<pat::Particle> fitBottom_, fitLight_, fitLepton_, fitNeutrino_;
00155 
00156   //-------------------------------------------
00157   // miscellaneous
00158   //-------------------------------------------
00159   std::string decay_;
00160   int jetCorrScheme_;
00161   double chi2Prob_;
00162   std::vector<double> scanValues_;
00163   double pTrueCombExist_, pTrueBJetSel_, pTrueBhadrSel_, pTrueJetComb_;
00164   double signalPur_, signalLRTot_;
00165   double sumDeltaRjp_, deltaRB_, deltaRL_;
00166   int changeBL_;
00167   bool bestSol_;
00168   //double jetMatchPur_;
00169 };
00170 
00171 #endif