File indexing completed on 2024-04-06 12:31:20
0001 #include "TopQuarkAnalysis/TopJetCombination/interface/TtFullLepHypothesis.h"
0002 #include "DataFormats/PatCandidates/interface/Particle.h"
0003
0004 class TtFullLepHypKinSolution : public TtFullLepHypothesis {
0005 public:
0006 explicit TtFullLepHypKinSolution(const edm::ParameterSet&);
0007
0008 private:
0009
0010 void buildKey() override { key_ = TtEvent::kKinSolution; };
0011
0012 void buildHypo(edm::Event& evt,
0013 const edm::Handle<std::vector<pat::Electron> >& elecs,
0014 const edm::Handle<std::vector<pat::Muon> >& mus,
0015 const edm::Handle<std::vector<pat::Jet> >& jets,
0016 const edm::Handle<std::vector<pat::MET> >& mets,
0017 std::vector<int>& match,
0018 const unsigned int iComb) override;
0019
0020
0021 edm::EDGetTokenT<std::vector<reco::LeafCandidate> > nusToken_;
0022 edm::EDGetTokenT<std::vector<reco::LeafCandidate> > nuBarsToken_;
0023 edm::EDGetTokenT<std::vector<double> > solWeightToken_;
0024 };
0025
0026 TtFullLepHypKinSolution::TtFullLepHypKinSolution(const edm::ParameterSet& cfg)
0027 : TtFullLepHypothesis(cfg),
0028 nusToken_(consumes<std::vector<reco::LeafCandidate> >(cfg.getParameter<edm::InputTag>("Neutrinos"))),
0029 nuBarsToken_(consumes<std::vector<reco::LeafCandidate> >(cfg.getParameter<edm::InputTag>("NeutrinoBars"))),
0030 solWeightToken_(consumes<std::vector<double> >(cfg.getParameter<edm::InputTag>("solutionWeight"))) {}
0031
0032 void TtFullLepHypKinSolution::buildHypo(edm::Event& evt,
0033 const edm::Handle<std::vector<pat::Electron> >& elecs,
0034 const edm::Handle<std::vector<pat::Muon> >& mus,
0035 const edm::Handle<std::vector<pat::Jet> >& jets,
0036 const edm::Handle<std::vector<pat::MET> >& mets,
0037 std::vector<int>& match,
0038 const unsigned int iComb) {
0039 edm::Handle<std::vector<double> > solWeight;
0040
0041 edm::Handle<std::vector<reco::LeafCandidate> > nus;
0042 edm::Handle<std::vector<reco::LeafCandidate> > nuBars;
0043
0044 evt.getByToken(solWeightToken_, solWeight);
0045
0046 evt.getByToken(nusToken_, nus);
0047 evt.getByToken(nuBarsToken_, nuBars);
0048
0049 if ((*solWeight)[iComb] < 0) {
0050
0051 return;
0052 }
0053
0054
0055
0056
0057 if (!jets->empty()) {
0058 b_ = makeCandidate(jets, match[0], jetCorrectionLevel_);
0059 bBar_ = makeCandidate(jets, match[1], jetCorrectionLevel_);
0060 }
0061
0062
0063
0064 if (!elecs->empty() && match[2] >= 0)
0065 leptonBar_ = makeCandidate(elecs, match[2]);
0066
0067 if (!elecs->empty() && match[3] >= 0)
0068 lepton_ = makeCandidate(elecs, match[3]);
0069
0070 if (!mus->empty() && match[4] >= 0 && match[2] < 0)
0071 leptonBar_ = makeCandidate(mus, match[4]);
0072
0073
0074
0075
0076 else if (!mus->empty() && match[4] >= 0)
0077 lepton_ = makeCandidate(mus, match[4]);
0078
0079 if (!mus->empty() && match[5] >= 0 && match[3] < 0)
0080 lepton_ = makeCandidate(mus, match[5]);
0081
0082
0083
0084
0085 else if (!mus->empty() && match[5] >= 0)
0086 leptonBar_ = makeCandidate(mus, match[5]);
0087
0088
0089
0090
0091 if (!nus->empty())
0092 neutrino_ = makeCandidate(nus, iComb);
0093
0094 if (!nuBars->empty())
0095 neutrinoBar_ = makeCandidate(nuBars, iComb);
0096 }
0097
0098 #include "FWCore/Framework/interface/MakerMacros.h"
0099 DEFINE_FWK_MODULE(TtFullLepHypKinSolution);