CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/AnalysisDataFormats/TopObjects/interface/TtHadEvtSolution.h

Go to the documentation of this file.
00001 #ifndef TopObjects_TtHadEvtSolution_h
00002 #define TopObjects_TtHadEvtSolution_h
00003 //
00004 // $Id: TtHadEvtSolution.h,v 1.13 2013/04/19 22:13:23 wmtan Exp $
00005 // adapted TtSemiEvtSolution.h,v 1.14 2007/07/06 03:07:47 lowette Exp 
00006 // for fully hadronic channel
00007 
00008 #include <vector>
00009 #include <string>
00010 
00011 #include "DataFormats/Common/interface/Handle.h"
00012 #include "DataFormats/Common/interface/RefProd.h"
00013 #include "DataFormats/Common/interface/Ref.h"
00014 
00015 #include "DataFormats/Candidate/interface/Particle.h"
00016 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00017 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
00018 
00019 #include "DataFormats/PatCandidates/interface/Particle.h"
00020 #include "DataFormats/PatCandidates/interface/Jet.h"
00021 
00022 class TtHadEvtSolution {
00023   
00024   friend class TtHadEvtSolutionMaker;
00025   friend class TtFullHadKinFitter;
00026   friend class TtHadLRJetCombObservables;
00027   friend class TtHadLRJetCombCalc;
00028   /*
00029     friend class TtHadLRSignalSelObservables;
00030     friend class TtHadLRSignalSelCalc;
00031   */
00032   
00033   public:
00034   
00035   TtHadEvtSolution();
00036   virtual ~TtHadEvtSolution();     
00037   
00038   //-------------------------------------------
00039   // get calibrated base objects 
00040   //-------------------------------------------
00041   pat::Jet getHadb() const;
00042   pat::Jet getHadp() const;
00043   pat::Jet getHadq() const;
00044   pat::Jet getHadbbar() const;
00045   pat::Jet getHadj() const;
00046   pat::Jet getHadk() const;
00047   
00048   //-------------------------------------------
00049   // get the matched gen particles
00050   //-------------------------------------------
00051   const edm::RefProd<TtGenEvent> & getGenEvent() const { return theGenEvt_; };
00052   const reco::GenParticle * getGenHadb() const { if (!theGenEvt_) return 0; else return theGenEvt_->b(); };
00053   const reco::GenParticle * getGenHadbbar() const { if (!theGenEvt_) return 0; else return theGenEvt_->bBar(); };
00054   const reco::GenParticle * getGenHadp() const { if (!theGenEvt_) return 0; else return theGenEvt_->daughterQuarkOfWPlus(); };
00055   const reco::GenParticle * getGenHadq() const { if (!theGenEvt_) return 0; else return theGenEvt_->daughterQuarkBarOfWPlus(); };
00056   const reco::GenParticle * getGenHadj() const { if (!theGenEvt_) return 0; else return theGenEvt_->daughterQuarkOfWMinus(); };
00057   const reco::GenParticle * getGenHadk() const { if (!theGenEvt_) return 0; else return theGenEvt_->daughterQuarkBarOfWMinus(); };
00058   
00059   //-------------------------------------------
00060   // get (un-)/calibrated reco objects
00061   //-------------------------------------------
00062   reco::Particle getRecHadt() const;
00063   reco::Particle getRecHadtbar() const;
00064   reco::Particle getRecHadW_plus() const;     
00065   reco::Particle getRecHadW_minus() const;       
00066   
00067   pat::Jet getRecHadb() const { return this->getHadb().correctedJet("RAW"); };
00068   pat::Jet getRecHadbbar() const { return this->getHadbbar().correctedJet("RAW"); };
00069   pat::Jet getRecHadp() const { return this->getHadp().correctedJet("RAW"); };
00070   pat::Jet getRecHadq() const { return this->getHadq().correctedJet("RAW"); };
00071   pat::Jet getRecHadj() const { return this->getHadj().correctedJet("RAW"); };
00072   pat::Jet getRecHadk() const { return this->getHadk().correctedJet("RAW"); };
00073   
00074   reco::Particle getCalHadt() const;
00075   reco::Particle getCalHadtbar() const;
00076   reco::Particle getCalHadW_plus() const;
00077   reco::Particle getCalHadW_minus() const;
00078   pat::Jet getCalHadb() const { return this->getHadb(); };
00079   pat::Jet getCalHadbbar() const { return this->getHadbbar(); };
00080   pat::Jet getCalHadp() const { return this->getHadp(); };
00081   pat::Jet getCalHadq() const { return this->getHadq(); };
00082   pat::Jet getCalHadj() const { return this->getHadj(); };
00083   pat::Jet getCalHadk() const { return this->getHadk(); };
00084 
00085   //-------------------------------------------
00086   // get objects from kinematic fit
00087   //-------------------------------------------  
00088   reco::Particle getFitHadt() const;
00089   reco::Particle getFitHadtbar() const;
00090   reco::Particle getFitHadW_plus() const;
00091   reco::Particle getFitHadW_minus() const;
00092   pat::Particle getFitHadb() const { return (fitHadb_.size()>0 ? fitHadb_.front() : pat::Particle()); };
00093   pat::Particle getFitHadbbar() const { return (fitHadbbar_.size()>0 ? fitHadbbar_.front() : pat::Particle()); };
00094   pat::Particle getFitHadp() const { return (fitHadp_.size()>0 ? fitHadp_.front() : pat::Particle()); };
00095   pat::Particle getFitHadq() const { return (fitHadq_.size()>0 ? fitHadq_.front() : pat::Particle()); };
00096   pat::Particle getFitHadj() const { return (fitHadj_.size()>0 ? fitHadj_.front() : pat::Particle()); };
00097   pat::Particle getFitHadk() const { return (fitHadk_.size()>0 ? fitHadk_.front() : pat::Particle()); };
00098 
00099   //-------------------------------------------  
00100   // get the selected hadronic decay chain 
00101   //-------------------------------------------
00102   std::string getDecay() const { return decay_; }
00103 
00104   //-------------------------------------------  
00105   // get info on the matching
00106   //-------------------------------------------  
00107   double getMCBestSumAngles() const { return sumAnglejp_; };
00108   double getMCBestAngleHadp() const { return angleHadp_; };
00109   double getMCBestAngleHadq() const { return angleHadq_; };
00110   double getMCBestAngleHadj() const { return angleHadj_; };
00111   double getMCBestAngleHadk() const { return angleHadk_; };
00112   double getMCBestAngleHadb() const { return angleHadb_; };
00113   double getMCBestAngleHadbbar() const { return angleHadbbar_; };
00114   int getMCChangeW1Q() const { return changeW1Q_; };     
00115   int getMCChangeW2Q() const { return changeW2Q_;}; 
00116 
00117   //-------------------------------------------  
00118   // get selected kinfit parametrisations of 
00119   //each type of object 
00120   //-------------------------------------------  
00121   int getJetParametrisation() const { return jetParam_; }
00122 
00123   //-------------------------------------------  
00124   // get the prob of the chi2 value resulting 
00125   // from the kinematic fit added chi2 for all 
00126   // fits
00127   //-------------------------------------------  
00128   double getProbChi2() const { return probChi2_; }
00129 
00130   //-------------------------------------------  
00131   // get info on the outcome of the signal 
00132   // selection LR
00133   //-------------------------------------------  
00134   double getLRSignalEvtObsVal(unsigned int) const;
00135   double getLRSignalEvtLRval() const { return lrSignalEvtLRval_; }
00136   double getLRSignalEvtProb() const { return lrSignalEvtProb_; }
00137 
00138   //-------------------------------------------  
00139   // get info on the outcome of the different 
00140   //jet combination methods
00141   //-------------------------------------------  
00142   int getMCBestJetComb() const { return mcBestJetComb_; }
00143   int getSimpleBestJetComb() const { return simpleBestJetComb_; }
00144   int getLRBestJetComb() const { return lrBestJetComb_; }
00145   double getLRJetCombObsVal(unsigned int) const;
00146   double getLRJetCombLRval() const { return lrJetCombLRval_; }
00147   double getLRJetCombProb() const { return lrJetCombProb_; }
00148   
00149   //protected: seems to cause compile error, check!!!
00150   
00151   //-------------------------------------------  
00152   // set the generated event
00153   //-------------------------------------------  
00154   void setGenEvt(const edm::Handle<TtGenEvent> & aGenEvt);
00155 
00156   //-------------------------------------------
00157   // set the basic objects
00158   //-------------------------------------------  
00159   void setJetCorrectionScheme(int scheme) { jetCorrScheme_ = scheme; };
00160   void setHadp(const edm::Handle<std::vector<pat::Jet> > & jet, int i)
00161   { hadp_ = edm::Ref<std::vector<pat::Jet> >(jet, i); }
00162   void setHadq(const edm::Handle<std::vector<pat::Jet> > & jet, int i)
00163   { hadq_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00164   void setHadj(const edm::Handle<std::vector<pat::Jet> > & jet, int i)
00165   { hadj_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00166   void setHadk(const edm::Handle<std::vector<pat::Jet> > & jet, int i)
00167   { hadk_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00168   void setHadb(const edm::Handle<std::vector<pat::Jet> > & jet, int i)
00169   { hadb_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00170   void setHadbbar(const edm::Handle<std::vector<pat::Jet> > & jet, int i)
00171   { hadbbar_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00172 
00173   //-------------------------------------------
00174   // set the fitted objects 
00175   //-------------------------------------------
00176   void setFitHadp(const pat::Particle & aFitHadp) { fitHadp_.clear(); fitHadp_.push_back(aFitHadp); };
00177   void setFitHadq(const pat::Particle & aFitHadq) { fitHadq_.clear(); fitHadq_.push_back(aFitHadq); };
00178   void setFitHadj(const pat::Particle & aFitHadj) { fitHadj_.clear(); fitHadj_.push_back(aFitHadj); };
00179   void setFitHadk(const pat::Particle & aFitHadk) { fitHadk_.clear(); fitHadk_.push_back(aFitHadk); };
00180   void setFitHadb(const pat::Particle & aFitHadb) { fitHadb_.clear(); fitHadb_.push_back(aFitHadb); };
00181   void setFitHadbbar(const pat::Particle & aFitHadbbar) { fitHadbbar_.clear(); fitHadbbar_.push_back(aFitHadbbar); };
00182 
00183   //-------------------------------------------
00184   // set matching info
00185   //-------------------------------------------
00186   void setMCBestSumAngles(double sdr) { sumAnglejp_ = sdr; };
00187   void setMCBestAngleHadp(double adr) { angleHadp_ = adr; };
00188   void setMCBestAngleHadq(double adr) { angleHadq_ = adr; };
00189   void setMCBestAngleHadj(double adr) { angleHadj_ = adr; };
00190   void setMCBestAngleHadk(double adr) { angleHadk_ = adr; };
00191   void setMCBestAngleHadb(double adr) { angleHadb_ = adr; };
00192   void setMCBestAngleHadbbar(double adr) { angleHadbbar_ = adr; };
00193   void setMCChangeW1Q(int w1q) { changeW1Q_ = w1q; };
00194   void setMCChangeW2Q(int w2q) { changeW2Q_ = w2q; };
00195 
00196   //-------------------------------------------
00197   // methods to set the kinfit parametrisations 
00198   //of each type of object 
00199   //-------------------------------------------
00200   void setJetParametrisation(int jp) { jetParam_ = jp; };
00201 
00202   //-------------------------------------------
00203   // method to set the prob. of the chi2 value 
00204   //resulting from the kinematic fit 
00205   //-------------------------------------------
00206   void setProbChi2(double c) { probChi2_ = c; };
00207 
00208   //-------------------------------------------
00209   // methods to set the outcome of the different 
00210   // jet combination methods
00211   //-------------------------------------------
00212   void setMCBestJetComb(int mcbs) { mcBestJetComb_ = mcbs; };
00213   void setSimpleBestJetComb(int sbs) { simpleBestJetComb_ = sbs; };
00214   void setLRBestJetComb(int lrbs) { lrBestJetComb_ = lrbs; };
00215   void setLRJetCombObservables(const std::vector<std::pair<unsigned int, double> >& varval);
00216   void setLRJetCombLRval(double clr) {lrJetCombLRval_ = clr;};
00217   void setLRJetCombProb(double plr) {lrJetCombProb_ = plr;};
00218 
00219   //-------------------------------------------
00220   // methods to set the outcome of the signal 
00221   // selection LR
00222   //-------------------------------------------
00223   void setLRSignalEvtObservables(const std::vector<std::pair<unsigned int, double> >& varval);
00224   void setLRSignalEvtLRval(double clr) {lrSignalEvtLRval_ = clr;};
00225   void setLRSignalEvtProb(double plr) {lrSignalEvtProb_ = plr;};
00226   
00227  private:
00228 
00229   //-------------------------------------------  
00230   // particle content
00231   //-------------------------------------------  
00232   edm::RefProd<TtGenEvent> theGenEvt_;
00233   edm::Ref<std::vector<pat::Jet> > hadb_, hadp_, hadq_, hadbbar_,hadj_, hadk_;
00234   std::vector<pat::Particle> fitHadb_, fitHadp_, fitHadq_, fitHadbbar_, fitHadj_, fitHadk_;
00235   
00236   std::string decay_;
00237   int jetCorrScheme_;
00238   double sumAnglejp_, angleHadp_, angleHadq_, angleHadb_, angleHadbbar_, angleHadj_ , angleHadk_;
00239   int changeW1Q_, changeW2Q_;
00240   int jetParam_;
00241   double probChi2_;
00242   int mcBestJetComb_, simpleBestJetComb_, lrBestJetComb_;
00243   double lrJetCombLRval_, lrJetCombProb_;
00244   double lrSignalEvtLRval_, lrSignalEvtProb_;
00245   std::vector<std::pair<unsigned int, double> > lrJetCombVarVal_;
00246   std::vector<std::pair<unsigned int, double> > lrSignalEvtVarVal_;
00247 };
00248 
00249 #endif