File indexing completed on 2024-04-06 11:57:31
0001
0002
0003
0004 #ifndef TopObjects_TtSemiEvtSolution_h
0005 #define TopObjects_TtSemiEvtSolution_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 #include "DataFormats/Candidate/interface/ShallowClonePtrCandidate.h"
0026 #include "DataFormats/Candidate/interface/CompositeCandidate.h"
0027
0028
0029
0030
0031 class TtSemiEvtSolution {
0032 friend class TtSemiEvtSolutionMaker;
0033 friend class TtSemiLepKinFitter;
0034 friend class TtSemiLepHitFit;
0035 friend class TtSemiLRSignalSelObservables;
0036 friend class TtSemiLRSignalSelCalc;
0037 friend class TtSemiLRJetCombObservables;
0038 friend class TtSemiLRJetCombCalc;
0039
0040 public:
0041 TtSemiEvtSolution();
0042 virtual ~TtSemiEvtSolution();
0043
0044
0045
0046
0047 pat::Jet getHadb() const;
0048 pat::Jet getHadp() const;
0049 pat::Jet getHadq() const;
0050 pat::Jet getLepb() const;
0051 pat::Muon getMuon() const { return *muon_; };
0052 pat::Electron getElectron() const { return *electron_; };
0053 pat::MET getNeutrino() const { return *neutrino_; };
0054
0055
0056
0057
0058 const edm::RefProd<TtGenEvent>& getGenEvent() const { return theGenEvt_; };
0059 const reco::GenParticle* getGenHadt() const {
0060 if (!theGenEvt_)
0061 return nullptr;
0062 else
0063 return this->getGenEvent()->hadronicDecayTop();
0064 };
0065 const reco::GenParticle* getGenHadW() const {
0066 if (!theGenEvt_)
0067 return nullptr;
0068 else
0069 return this->getGenEvent()->hadronicDecayW();
0070 };
0071 const reco::GenParticle* getGenHadb() const {
0072 if (!theGenEvt_)
0073 return nullptr;
0074 else
0075 return this->getGenEvent()->hadronicDecayB();
0076 };
0077 const reco::GenParticle* getGenHadp() const {
0078 if (!theGenEvt_)
0079 return nullptr;
0080 else
0081 return this->getGenEvent()->hadronicDecayQuark();
0082 };
0083 const reco::GenParticle* getGenHadq() const {
0084 if (!theGenEvt_)
0085 return nullptr;
0086 else
0087 return this->getGenEvent()->hadronicDecayQuarkBar();
0088 };
0089 const reco::GenParticle* getGenLept() const {
0090 if (!theGenEvt_)
0091 return nullptr;
0092 else
0093 return this->getGenEvent()->leptonicDecayTop();
0094 };
0095 const reco::GenParticle* getGenLepW() const {
0096 if (!theGenEvt_)
0097 return nullptr;
0098 else
0099 return this->getGenEvent()->leptonicDecayW();
0100 };
0101 const reco::GenParticle* getGenLepb() const {
0102 if (!theGenEvt_)
0103 return nullptr;
0104 else
0105 return this->getGenEvent()->leptonicDecayB();
0106 };
0107 const reco::GenParticle* getGenLepl() const {
0108 if (!theGenEvt_)
0109 return nullptr;
0110 else
0111 return this->getGenEvent()->singleLepton();
0112 };
0113 const reco::GenParticle* getGenLepn() const {
0114 if (!theGenEvt_)
0115 return nullptr;
0116 else
0117 return this->getGenEvent()->singleNeutrino();
0118 };
0119
0120
0121
0122
0123 reco::Particle getRecHadt() const;
0124 reco::Particle getRecHadW() const;
0125 pat::Jet getRecHadb() const { return this->getHadb().correctedJet("RAW"); };
0126 pat::Jet getRecHadp() const { return this->getHadp().correctedJet("RAW"); };
0127 pat::Jet getRecHadq() const { return this->getHadq().correctedJet("RAW"); };
0128 reco::Particle getRecLept() const;
0129 reco::Particle getRecLepW() const;
0130 pat::Jet getRecLepb() const { return this->getLepb().correctedJet("RAW"); };
0131 pat::Muon getRecLepm() const { return this->getMuon(); };
0132 pat::Electron getRecLepe() const { return this->getElectron(); };
0133 pat::MET getRecLepn() const { return this->getNeutrino(); };
0134
0135
0136 reco::Particle getCalHadt() const;
0137 reco::Particle getCalHadW() const;
0138 pat::Jet getCalHadb() const { return this->getHadb(); };
0139 pat::Jet getCalHadp() const { return this->getHadp(); };
0140 pat::Jet getCalHadq() const { return this->getHadq(); };
0141 reco::Particle getCalLept() const;
0142 reco::Particle getCalLepW() const;
0143 pat::Jet getCalLepb() const { return this->getLepb(); };
0144 pat::Muon getCalLepm() const { return this->getMuon(); };
0145 pat::Electron getCalLepe() const { return this->getElectron(); };
0146 pat::MET getCalLepn() const { return this->getNeutrino(); };
0147
0148
0149
0150
0151 reco::Particle getFitHadt() const;
0152 reco::Particle getFitHadW() const;
0153 pat::Particle getFitHadb() const { return (!fitHadb_.empty() ? fitHadb_.front() : pat::Particle()); };
0154 pat::Particle getFitHadp() const { return (!fitHadp_.empty() ? fitHadp_.front() : pat::Particle()); };
0155 pat::Particle getFitHadq() const { return (!fitHadq_.empty() ? fitHadq_.front() : pat::Particle()); };
0156 reco::Particle getFitLept() const;
0157 reco::Particle getFitLepW() const;
0158 pat::Particle getFitLepb() const { return (!fitLepb_.empty() ? fitLepb_.front() : pat::Particle()); };
0159 pat::Particle getFitLepl() const { return (!fitLepl_.empty() ? fitLepl_.front() : pat::Particle()); };
0160 pat::Particle getFitLepn() const { return (!fitLepn_.empty() ? fitLepn_.front() : pat::Particle()); };
0161
0162
0163
0164
0165 std::string getDecay() const { return decay_; }
0166
0167
0168
0169
0170 double getMCBestSumAngles() const { return sumAnglejp_; };
0171 double getMCBestAngleHadp() const { return angleHadp_; };
0172 double getMCBestAngleHadq() const { return angleHadq_; };
0173 double getMCBestAngleHadb() const { return angleHadb_; };
0174 double getMCBestAngleLepb() const { return angleLepb_; };
0175 int getMCChangeWQ() const { return changeWQ_; };
0176
0177
0178
0179
0180
0181 int getJetParametrisation() const { return jetParam_; }
0182 int getLeptonParametrisation() const { return leptonParam_; }
0183 int getNeutrinoParametrisation() const { return neutrinoParam_; }
0184
0185
0186
0187
0188
0189 double getProbChi2() const { return probChi2_; }
0190
0191
0192
0193
0194
0195 double getLRSignalEvtObsVal(unsigned int) const;
0196 double getLRSignalEvtLRval() const { return lrSignalEvtLRval_; }
0197 double getLRSignalEvtProb() const { return lrSignalEvtProb_; }
0198
0199
0200
0201
0202
0203 int getMCBestJetComb() const { return mcBestJetComb_; }
0204 int getSimpleBestJetComb() const { return simpleBestJetComb_; }
0205 int getLRBestJetComb() const { return lrBestJetComb_; }
0206 double getLRJetCombObsVal(unsigned int) const;
0207 double getLRJetCombLRval() const { return lrJetCombLRval_; }
0208 double getLRJetCombProb() const { return lrJetCombProb_; }
0209
0210
0211
0212
0213 const reco::CompositeCandidate& getRecoHyp() const { return recoHyp_; }
0214 const reco::CompositeCandidate& getFitHyp() const { return fitHyp_; }
0215 const reco::CompositeCandidate& getMCHyp() const { return mcHyp_; }
0216
0217 protected:
0218
0219
0220
0221 void setGenEvt(const edm::Handle<TtGenEvent>& aGenEvt);
0222
0223
0224
0225
0226 void setJetCorrectionScheme(int scheme) { jetCorrScheme_ = scheme; };
0227 void setHadp(const edm::Handle<std::vector<pat::Jet> >& jet, int i) { hadp_ = edm::Ptr<pat::Jet>(jet, i); };
0228 void setHadq(const edm::Handle<std::vector<pat::Jet> >& jet, int i) { hadq_ = edm::Ptr<pat::Jet>(jet, i); };
0229 void setHadb(const edm::Handle<std::vector<pat::Jet> >& jet, int i) { hadb_ = edm::Ptr<pat::Jet>(jet, i); };
0230 void setLepb(const edm::Handle<std::vector<pat::Jet> >& jet, int i) { lepb_ = edm::Ptr<pat::Jet>(jet, i); };
0231 void setMuon(const edm::Handle<std::vector<pat::Muon> >& muon, int i) {
0232 muon_ = edm::Ptr<pat::Muon>(muon, i);
0233 decay_ = "muon";
0234 };
0235 void setElectron(const edm::Handle<std::vector<pat::Electron> >& elec, int i) {
0236 electron_ = edm::Ptr<pat::Electron>(elec, i);
0237 decay_ = "electron";
0238 };
0239 void setNeutrino(const edm::Handle<std::vector<pat::MET> >& met, int i) { neutrino_ = edm::Ptr<pat::MET>(met, i); };
0240
0241
0242
0243
0244 void setFitHadb(const pat::Particle& aFitHadb) {
0245 fitHadb_.clear();
0246 fitHadb_.push_back(aFitHadb);
0247 };
0248 void setFitHadp(const pat::Particle& aFitHadp) {
0249 fitHadp_.clear();
0250 fitHadp_.push_back(aFitHadp);
0251 };
0252 void setFitHadq(const pat::Particle& aFitHadq) {
0253 fitHadq_.clear();
0254 fitHadq_.push_back(aFitHadq);
0255 };
0256 void setFitLepb(const pat::Particle& aFitLepb) {
0257 fitLepb_.clear();
0258 fitLepb_.push_back(aFitLepb);
0259 };
0260 void setFitLepl(const pat::Particle& aFitLepl) {
0261 fitLepl_.clear();
0262 fitLepl_.push_back(aFitLepl);
0263 };
0264 void setFitLepn(const pat::Particle& aFitLepn) {
0265 fitLepn_.clear();
0266 fitLepn_.push_back(aFitLepn);
0267 };
0268
0269
0270
0271
0272 void setMCBestSumAngles(double sdr) { sumAnglejp_ = sdr; };
0273 void setMCBestAngleHadp(double adr) { angleHadp_ = adr; };
0274 void setMCBestAngleHadq(double adr) { angleHadq_ = adr; };
0275 void setMCBestAngleHadb(double adr) { angleHadb_ = adr; };
0276 void setMCBestAngleLepb(double adr) { angleLepb_ = adr; };
0277 void setMCChangeWQ(int wq) { changeWQ_ = wq; };
0278
0279
0280
0281
0282
0283 void setJetParametrisation(int jp) { jetParam_ = jp; };
0284 void setLeptonParametrisation(int lp) { leptonParam_ = lp; };
0285 void setNeutrinoParametrisation(int mp) { neutrinoParam_ = mp; };
0286
0287
0288
0289
0290
0291 void setProbChi2(double c) { probChi2_ = c; };
0292
0293
0294
0295
0296
0297 void setMCBestJetComb(int mcbs) { mcBestJetComb_ = mcbs; };
0298 void setSimpleBestJetComb(int sbs) { simpleBestJetComb_ = sbs; };
0299 void setLRBestJetComb(int lrbs) { lrBestJetComb_ = lrbs; };
0300 void setLRJetCombObservables(const std::vector<std::pair<unsigned int, double> >& varval);
0301 void setLRJetCombLRval(double clr) { lrJetCombLRval_ = clr; };
0302 void setLRJetCombProb(double plr) { lrJetCombProb_ = plr; };
0303
0304
0305
0306
0307 void setLRSignalEvtObservables(const std::vector<std::pair<unsigned int, double> >& varval);
0308 void setLRSignalEvtLRval(double clr) { lrSignalEvtLRval_ = clr; };
0309 void setLRSignalEvtProb(double plr) { lrSignalEvtProb_ = plr; };
0310
0311 private:
0312
0313
0314
0315 edm::RefProd<TtGenEvent> theGenEvt_;
0316 edm::Ptr<pat::Jet> hadb_, hadp_, hadq_, lepb_;
0317 edm::Ptr<pat::Muon> muon_;
0318 edm::Ptr<pat::Electron> electron_;
0319 edm::Ptr<pat::MET> neutrino_;
0320 std::vector<pat::Particle> fitHadb_, fitHadp_, fitHadq_;
0321 std::vector<pat::Particle> fitLepb_, fitLepl_, fitLepn_;
0322
0323 reco::CompositeCandidate mcHyp_;
0324 reco::CompositeCandidate recoHyp_;
0325 reco::CompositeCandidate fitHyp_;
0326
0327 void setupHyp();
0328
0329 std::string decay_;
0330 int jetCorrScheme_;
0331 double sumAnglejp_, angleHadp_, angleHadq_, angleHadb_, angleLepb_;
0332 int changeWQ_;
0333 int jetParam_, leptonParam_, neutrinoParam_;
0334 double probChi2_;
0335 int mcBestJetComb_, simpleBestJetComb_, lrBestJetComb_;
0336 double lrJetCombLRval_, lrJetCombProb_;
0337 double lrSignalEvtLRval_, lrSignalEvtProb_;
0338 std::vector<std::pair<unsigned int, double> > lrJetCombVarVal_;
0339 std::vector<std::pair<unsigned int, double> > lrSignalEvtVarVal_;
0340 };
0341
0342 #endif