00001
00002
00003
00004
00005 #ifndef TopObjects_TtDilepEvtSolution_h
00006 #define TopObjects_TtDilepEvtSolution_h
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/Electron.h"
00021 #include "DataFormats/PatCandidates/interface/Muon.h"
00022 #include "DataFormats/PatCandidates/interface/Tau.h"
00023 #include "DataFormats/PatCandidates/interface/Jet.h"
00024 #include "DataFormats/PatCandidates/interface/MET.h"
00025
00026 class TtDilepEvtSolution {
00027
00028 friend class TtFullLepKinSolver;
00029 friend class TtDilepEvtSolutionMaker;
00030 friend class TtDilepLRSignalSelObservables;
00031 friend class TtLRSignalSelCalc;
00032
00033 public:
00034
00035 TtDilepEvtSolution();
00036 virtual ~TtDilepEvtSolution();
00037
00038
00039
00040
00041 pat::Jet getJetB() const;
00042 pat::Jet getJetBbar() const;
00043 pat::Electron getElectronp() const { return *elecp_; };
00044 pat::Electron getElectronm() const { return *elecm_; };
00045 pat::Muon getMuonp() const { return *muonp_; };
00046 pat::Muon getMuonm() const { return *muonm_; };
00047 pat::Tau getTaup() const { return *taup_; };
00048 pat::Tau getTaum() const { return *taum_; };
00049 pat::MET getMET() const { return *met_; };
00050
00051
00052
00053
00054 const edm::RefProd<TtGenEvent> & getGenEvent() const { return theGenEvt_; };
00055 const reco::GenParticle * getGenT() const { if (!theGenEvt_) return 0; else return theGenEvt_->top(); };
00056 const reco::GenParticle * getGenWp() const { if (!theGenEvt_) return 0; else return theGenEvt_->wPlus(); };
00057 const reco::GenParticle * getGenB() const { if (!theGenEvt_) return 0; else return theGenEvt_->b(); };
00058 const reco::GenParticle * getGenLepp() const { if (!theGenEvt_) return 0; else return theGenEvt_->leptonBar(); };
00059 const reco::GenParticle * getGenN() const { if (!theGenEvt_) return 0; else return theGenEvt_->neutrino(); };
00060 const reco::GenParticle * getGenTbar() const { if (!theGenEvt_) return 0; else return theGenEvt_->topBar(); };
00061 const reco::GenParticle * getGenWm() const { if (!theGenEvt_) return 0; else return theGenEvt_->wMinus(); };
00062 const reco::GenParticle * getGenBbar() const { if (!theGenEvt_) return 0; else return theGenEvt_->bBar(); };
00063 const reco::GenParticle * getGenLepm() const { if (!theGenEvt_) return 0; else return theGenEvt_->lepton(); };
00064 const reco::GenParticle * getGenNbar() const { if (!theGenEvt_) return 0; else return theGenEvt_->neutrinoBar(); };
00065
00066
00067
00068
00069 pat::JetType getRecJetB() const { return this->getJetB().correctedJet("RAW"); };
00070 pat::Jet getCalJetB() const { return this->getJetB(); };
00071 pat::JetType getRecJetBbar() const { return this->getJetBbar().correctedJet("RAW"); };
00072 pat::Jet getCalJetBbar() const { return this->getJetBbar(); };
00073
00074
00075
00076
00077 std::string getWpDecay() const { return wpDecay_; }
00078 std::string getWmDecay() const { return wmDecay_; }
00079
00080
00081
00082
00083 double getJetResidual() const;
00084 double getLeptonResidual() const;
00085 double getFullResidual() const { return getJetResidual()+getFullResidual(); }
00086 bool getBestSol() const { return bestSol_; }
00087 double getRecTopMass() const {return topmass_; }
00088 double getRecWeightMax() const {return weightmax_; }
00089
00090
00091
00092
00093
00094 reco::Particle getLeptPos() const;
00095
00096
00097
00098
00099
00100 reco::Particle getLeptNeg() const;
00101
00102
00103
00104
00105
00106 double getLRSignalEvtObsVal(unsigned int) const;
00107 double getLRSignalEvtLRval() const { return lrSignalEvtLRval_; }
00108 double getLRSignalEvtProb() const { return lrSignalEvtProb_; }
00109
00110 protected:
00111
00112
00113
00114
00115 void setGenEvt(const edm::Handle<TtGenEvent>&);
00116
00117
00118
00119
00120 void setJetCorrectionScheme(int jetCorrScheme)
00121 { jetCorrScheme_ = jetCorrScheme; };
00122 void setB(const edm::Handle<std::vector<pat::Jet> >& jet, int i)
00123 { jetB_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00124 void setBbar(const edm::Handle<std::vector<pat::Jet> >& jet, int i)
00125 { jetBbar_ = edm::Ref<std::vector<pat::Jet> >(jet, i); };
00126 void setMuonp(const edm::Handle<std::vector<pat::Muon> >& muon, int i)
00127 { muonp_ = edm::Ref<std::vector<pat::Muon> >(muon, i); wpDecay_ = "muon"; };
00128 void setMuonm(const edm::Handle<std::vector<pat::Muon> >& muon, int i)
00129 { muonm_ = edm::Ref<std::vector<pat::Muon> >(muon, i); wmDecay_ = "muon"; }
00130 void setTaup(const edm::Handle<std::vector<pat::Tau> >& tau, int i)
00131 { taup_ = edm::Ref<std::vector<pat::Tau> >(tau, i); wpDecay_ = "tau"; }
00132 void setTaum(const edm::Handle<std::vector<pat::Tau> >& tau, int i)
00133 { taum_ = edm::Ref<std::vector<pat::Tau> >(tau, i); wmDecay_ = "tau"; }
00134 void setElectronp(const edm::Handle<std::vector<pat::Electron> >& elec, int i)
00135 { elecp_ = edm::Ref<std::vector<pat::Electron> >(elec, i); wpDecay_ = "electron"; };
00136 void setElectronm(const edm::Handle<std::vector<pat::Electron> >& elec, int i)
00137 { elecm_ = edm::Ref<std::vector<pat::Electron> >(elec, i); wmDecay_ = "electron"; };
00138 void setMET(const edm::Handle<std::vector<pat::MET> >& met, int i)
00139 { met_ = edm::Ref<std::vector<pat::MET> >(met, i); };
00140
00141
00142
00143
00144 void setBestSol(bool bs) { bestSol_ = bs; };
00145 void setRecTopMass(double mass) { topmass_ = mass; };
00146 void setRecWeightMax(double wgt) { weightmax_ = wgt; };
00147
00148
00149
00150
00151 void setLRSignalEvtObservables(std::vector<std::pair<unsigned int, double> >);
00152 void setLRSignalEvtLRval(double clr) {lrSignalEvtLRval_ = clr;};
00153 void setLRSignalEvtProb(double plr) {lrSignalEvtProb_ = plr;};
00154
00155 private:
00156
00157
00158
00159
00160 edm::RefProd<TtGenEvent> theGenEvt_;
00161 edm::Ref<std::vector<pat::Electron> > elecp_, elecm_;
00162 edm::Ref<std::vector<pat::Muon> > muonp_, muonm_;
00163 edm::Ref<std::vector<pat::Tau> > taup_, taum_;
00164 edm::Ref<std::vector<pat::Jet> > jetB_, jetBbar_;
00165 edm::Ref<std::vector<pat::MET> > met_;
00166
00167
00168
00169
00170 int jetCorrScheme_;
00171 std::string wpDecay_;
00172 std::string wmDecay_;
00173 bool bestSol_;
00174 double topmass_;
00175 double weightmax_;
00176
00177 double lrSignalEvtLRval_, lrSignalEvtProb_;
00178 std::vector<std::pair<unsigned int, double> > lrSignalEvtVarVal_;
00179 };
00180
00181 #endif