Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef TopObjects_TtHadEvtSolution_h
0002 #define TopObjects_TtHadEvtSolution_h
0003 //
0004 // adapted TtSemiEvtSolution.h,v 1.14 2007/07/06 03:07:47 lowette Exp
0005 // for fully hadronic channel
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/Jet.h"
0020 
0021 class TtHadEvtSolution {
0022   friend class TtHadEvtSolutionMaker;
0023   friend class TtFullHadKinFitter;
0024   friend class TtHadLRJetCombObservables;
0025   friend class TtHadLRJetCombCalc;
0026   /*
0027     friend class TtHadLRSignalSelObservables;
0028     friend class TtHadLRSignalSelCalc;
0029   */
0030 
0031 public:
0032   TtHadEvtSolution();
0033   virtual ~TtHadEvtSolution();
0034 
0035   //-------------------------------------------
0036   // get calibrated base objects
0037   //-------------------------------------------
0038   pat::Jet getHadb() const;
0039   pat::Jet getHadp() const;
0040   pat::Jet getHadq() const;
0041   pat::Jet getHadbbar() const;
0042   pat::Jet getHadj() const;
0043   pat::Jet getHadk() const;
0044 
0045   //-------------------------------------------
0046   // get the matched gen particles
0047   //-------------------------------------------
0048   const edm::RefProd<TtGenEvent>& getGenEvent() const { return theGenEvt_; };
0049   const reco::GenParticle* getGenHadb() const {
0050     if (!theGenEvt_)
0051       return nullptr;
0052     else
0053       return theGenEvt_->b();
0054   };
0055   const reco::GenParticle* getGenHadbbar() const {
0056     if (!theGenEvt_)
0057       return nullptr;
0058     else
0059       return theGenEvt_->bBar();
0060   };
0061   const reco::GenParticle* getGenHadp() const {
0062     if (!theGenEvt_)
0063       return nullptr;
0064     else
0065       return theGenEvt_->daughterQuarkOfWPlus();
0066   };
0067   const reco::GenParticle* getGenHadq() const {
0068     if (!theGenEvt_)
0069       return nullptr;
0070     else
0071       return theGenEvt_->daughterQuarkBarOfWPlus();
0072   };
0073   const reco::GenParticle* getGenHadj() const {
0074     if (!theGenEvt_)
0075       return nullptr;
0076     else
0077       return theGenEvt_->daughterQuarkOfWMinus();
0078   };
0079   const reco::GenParticle* getGenHadk() const {
0080     if (!theGenEvt_)
0081       return nullptr;
0082     else
0083       return theGenEvt_->daughterQuarkBarOfWMinus();
0084   };
0085 
0086   //-------------------------------------------
0087   // get (un-)/calibrated reco objects
0088   //-------------------------------------------
0089   reco::Particle getRecHadt() const;
0090   reco::Particle getRecHadtbar() const;
0091   reco::Particle getRecHadW_plus() const;
0092   reco::Particle getRecHadW_minus() const;
0093 
0094   pat::Jet getRecHadb() const { return this->getHadb().correctedJet("RAW"); };
0095   pat::Jet getRecHadbbar() const { return this->getHadbbar().correctedJet("RAW"); };
0096   pat::Jet getRecHadp() const { return this->getHadp().correctedJet("RAW"); };
0097   pat::Jet getRecHadq() const { return this->getHadq().correctedJet("RAW"); };
0098   pat::Jet getRecHadj() const { return this->getHadj().correctedJet("RAW"); };
0099   pat::Jet getRecHadk() const { return this->getHadk().correctedJet("RAW"); };
0100 
0101   reco::Particle getCalHadt() const;
0102   reco::Particle getCalHadtbar() const;
0103   reco::Particle getCalHadW_plus() const;
0104   reco::Particle getCalHadW_minus() const;
0105   pat::Jet getCalHadb() const { return this->getHadb(); };
0106   pat::Jet getCalHadbbar() const { return this->getHadbbar(); };
0107   pat::Jet getCalHadp() const { return this->getHadp(); };
0108   pat::Jet getCalHadq() const { return this->getHadq(); };
0109   pat::Jet getCalHadj() const { return this->getHadj(); };
0110   pat::Jet getCalHadk() const { return this->getHadk(); };
0111 
0112   //-------------------------------------------
0113   // get objects from kinematic fit
0114   //-------------------------------------------
0115   reco::Particle getFitHadt() const;
0116   reco::Particle getFitHadtbar() const;
0117   reco::Particle getFitHadW_plus() const;
0118   reco::Particle getFitHadW_minus() const;
0119   pat::Particle getFitHadb() const { return (!fitHadb_.empty() ? fitHadb_.front() : pat::Particle()); };
0120   pat::Particle getFitHadbbar() const { return (!fitHadbbar_.empty() ? fitHadbbar_.front() : pat::Particle()); };
0121   pat::Particle getFitHadp() const { return (!fitHadp_.empty() ? fitHadp_.front() : pat::Particle()); };
0122   pat::Particle getFitHadq() const { return (!fitHadq_.empty() ? fitHadq_.front() : pat::Particle()); };
0123   pat::Particle getFitHadj() const { return (!fitHadj_.empty() ? fitHadj_.front() : pat::Particle()); };
0124   pat::Particle getFitHadk() const { return (!fitHadk_.empty() ? fitHadk_.front() : pat::Particle()); };
0125 
0126   //-------------------------------------------
0127   // get the selected hadronic decay chain
0128   //-------------------------------------------
0129   std::string getDecay() const { return decay_; }
0130 
0131   //-------------------------------------------
0132   // get info on the matching
0133   //-------------------------------------------
0134   double getMCBestSumAngles() const { return sumAnglejp_; };
0135   double getMCBestAngleHadp() const { return angleHadp_; };
0136   double getMCBestAngleHadq() const { return angleHadq_; };
0137   double getMCBestAngleHadj() const { return angleHadj_; };
0138   double getMCBestAngleHadk() const { return angleHadk_; };
0139   double getMCBestAngleHadb() const { return angleHadb_; };
0140   double getMCBestAngleHadbbar() const { return angleHadbbar_; };
0141   int getMCChangeW1Q() const { return changeW1Q_; };
0142   int getMCChangeW2Q() const { return changeW2Q_; };
0143 
0144   //-------------------------------------------
0145   // get selected kinfit parametrisations of
0146   //each type of object
0147   //-------------------------------------------
0148   int getJetParametrisation() const { return jetParam_; }
0149 
0150   //-------------------------------------------
0151   // get the prob of the chi2 value resulting
0152   // from the kinematic fit added chi2 for all
0153   // fits
0154   //-------------------------------------------
0155   double getProbChi2() const { return probChi2_; }
0156 
0157   //-------------------------------------------
0158   // get info on the outcome of the signal
0159   // selection LR
0160   //-------------------------------------------
0161   double getLRSignalEvtObsVal(unsigned int) const;
0162   double getLRSignalEvtLRval() const { return lrSignalEvtLRval_; }
0163   double getLRSignalEvtProb() const { return lrSignalEvtProb_; }
0164 
0165   //-------------------------------------------
0166   // get info on the outcome of the different
0167   //jet combination methods
0168   //-------------------------------------------
0169   int getMCBestJetComb() const { return mcBestJetComb_; }
0170   int getSimpleBestJetComb() const { return simpleBestJetComb_; }
0171   int getLRBestJetComb() const { return lrBestJetComb_; }
0172   double getLRJetCombObsVal(unsigned int) const;
0173   double getLRJetCombLRval() const { return lrJetCombLRval_; }
0174   double getLRJetCombProb() const { return lrJetCombProb_; }
0175 
0176   //protected: seems to cause compile error, check!!!
0177 
0178   //-------------------------------------------
0179   // set the generated event
0180   //-------------------------------------------
0181   void setGenEvt(const edm::Handle<TtGenEvent>& aGenEvt);
0182 
0183   //-------------------------------------------
0184   // set the basic objects
0185   //-------------------------------------------
0186   void setJetCorrectionScheme(int scheme) { jetCorrScheme_ = scheme; };
0187   void setHadp(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0188     hadp_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0189   }
0190   void setHadq(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0191     hadq_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0192   };
0193   void setHadj(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0194     hadj_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0195   };
0196   void setHadk(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0197     hadk_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0198   };
0199   void setHadb(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0200     hadb_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0201   };
0202   void setHadbbar(const edm::Handle<std::vector<pat::Jet> >& jet, int i) {
0203     hadbbar_ = edm::Ref<std::vector<pat::Jet> >(jet, i);
0204   };
0205 
0206   //-------------------------------------------
0207   // set the fitted objects
0208   //-------------------------------------------
0209   void setFitHadp(const pat::Particle& aFitHadp) {
0210     fitHadp_.clear();
0211     fitHadp_.push_back(aFitHadp);
0212   };
0213   void setFitHadq(const pat::Particle& aFitHadq) {
0214     fitHadq_.clear();
0215     fitHadq_.push_back(aFitHadq);
0216   };
0217   void setFitHadj(const pat::Particle& aFitHadj) {
0218     fitHadj_.clear();
0219     fitHadj_.push_back(aFitHadj);
0220   };
0221   void setFitHadk(const pat::Particle& aFitHadk) {
0222     fitHadk_.clear();
0223     fitHadk_.push_back(aFitHadk);
0224   };
0225   void setFitHadb(const pat::Particle& aFitHadb) {
0226     fitHadb_.clear();
0227     fitHadb_.push_back(aFitHadb);
0228   };
0229   void setFitHadbbar(const pat::Particle& aFitHadbbar) {
0230     fitHadbbar_.clear();
0231     fitHadbbar_.push_back(aFitHadbbar);
0232   };
0233 
0234   //-------------------------------------------
0235   // set matching info
0236   //-------------------------------------------
0237   void setMCBestSumAngles(double sdr) { sumAnglejp_ = sdr; };
0238   void setMCBestAngleHadp(double adr) { angleHadp_ = adr; };
0239   void setMCBestAngleHadq(double adr) { angleHadq_ = adr; };
0240   void setMCBestAngleHadj(double adr) { angleHadj_ = adr; };
0241   void setMCBestAngleHadk(double adr) { angleHadk_ = adr; };
0242   void setMCBestAngleHadb(double adr) { angleHadb_ = adr; };
0243   void setMCBestAngleHadbbar(double adr) { angleHadbbar_ = adr; };
0244   void setMCChangeW1Q(int w1q) { changeW1Q_ = w1q; };
0245   void setMCChangeW2Q(int w2q) { changeW2Q_ = w2q; };
0246 
0247   //-------------------------------------------
0248   // methods to set the kinfit parametrisations
0249   //of each type of object
0250   //-------------------------------------------
0251   void setJetParametrisation(int jp) { jetParam_ = jp; };
0252 
0253   //-------------------------------------------
0254   // method to set the prob. of the chi2 value
0255   //resulting from the kinematic fit
0256   //-------------------------------------------
0257   void setProbChi2(double c) { probChi2_ = c; };
0258 
0259   //-------------------------------------------
0260   // methods to set the outcome of the different
0261   // jet combination methods
0262   //-------------------------------------------
0263   void setMCBestJetComb(int mcbs) { mcBestJetComb_ = mcbs; };
0264   void setSimpleBestJetComb(int sbs) { simpleBestJetComb_ = sbs; };
0265   void setLRBestJetComb(int lrbs) { lrBestJetComb_ = lrbs; };
0266   void setLRJetCombObservables(const std::vector<std::pair<unsigned int, double> >& varval);
0267   void setLRJetCombLRval(double clr) { lrJetCombLRval_ = clr; };
0268   void setLRJetCombProb(double plr) { lrJetCombProb_ = plr; };
0269 
0270   //-------------------------------------------
0271   // methods to set the outcome of the signal
0272   // selection LR
0273   //-------------------------------------------
0274   void setLRSignalEvtObservables(const std::vector<std::pair<unsigned int, double> >& varval);
0275   void setLRSignalEvtLRval(double clr) { lrSignalEvtLRval_ = clr; };
0276   void setLRSignalEvtProb(double plr) { lrSignalEvtProb_ = plr; };
0277 
0278 private:
0279   //-------------------------------------------
0280   // particle content
0281   //-------------------------------------------
0282   edm::RefProd<TtGenEvent> theGenEvt_;
0283   edm::Ref<std::vector<pat::Jet> > hadb_, hadp_, hadq_, hadbbar_, hadj_, hadk_;
0284   std::vector<pat::Particle> fitHadb_, fitHadp_, fitHadq_, fitHadbbar_, fitHadj_, fitHadk_;
0285 
0286   std::string decay_;
0287   int jetCorrScheme_;
0288   double sumAnglejp_, angleHadp_, angleHadq_, angleHadb_, angleHadbbar_, angleHadj_, angleHadk_;
0289   int changeW1Q_, changeW2Q_;
0290   int jetParam_;
0291   double probChi2_;
0292   int mcBestJetComb_, simpleBestJetComb_, lrBestJetComb_;
0293   double lrJetCombLRval_, lrJetCombProb_;
0294   double lrSignalEvtLRval_, lrSignalEvtProb_;
0295   std::vector<std::pair<unsigned int, double> > lrJetCombVarVal_;
0296   std::vector<std::pair<unsigned int, double> > lrSignalEvtVarVal_;
0297 };
0298 
0299 #endif