File indexing completed on 2024-04-06 12:31:20
0001 #include "TopQuarkAnalysis/TopJetCombination/interface/TtFullLepHypothesis.h"
0002
0003 #include "AnalysisDataFormats/TopObjects/interface/TtFullLepEvtPartons.h"
0004 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0005 #include "DataFormats/Math/interface/deltaR.h"
0006
0007 class TtFullLepHypGenMatch : public TtFullLepHypothesis {
0008 public:
0009 explicit TtFullLepHypGenMatch(const edm::ParameterSet&);
0010
0011 private:
0012
0013 void buildKey() override { key_ = TtFullLeptonicEvent::kGenMatch; };
0014
0015 void buildHypo(edm::Event& evt,
0016 const edm::Handle<std::vector<pat::Electron> >& elecs,
0017 const edm::Handle<std::vector<pat::Muon> >& mus,
0018 const edm::Handle<std::vector<pat::Jet> >& jets,
0019 const edm::Handle<std::vector<pat::MET> >& mets,
0020 std::vector<int>& match,
0021 const unsigned int iComb) override;
0022
0023 template <typename O>
0024 int findMatchingLepton(const reco::GenParticle*, const edm::Handle<std::vector<O> >&);
0025 void buildMatchingNeutrinos(edm::Event&, const edm::Handle<std::vector<pat::MET> >&);
0026
0027 protected:
0028 edm::EDGetTokenT<TtGenEvent> genEvtToken_;
0029 };
0030
0031 TtFullLepHypGenMatch::TtFullLepHypGenMatch(const edm::ParameterSet& cfg)
0032 : TtFullLepHypothesis(cfg), genEvtToken_(consumes<TtGenEvent>(edm::InputTag("genEvt"))) {}
0033
0034 void TtFullLepHypGenMatch::buildHypo(edm::Event& evt,
0035 const edm::Handle<std::vector<pat::Electron> >& elecs,
0036 const edm::Handle<std::vector<pat::Muon> >& mus,
0037 const edm::Handle<std::vector<pat::Jet> >& jets,
0038 const edm::Handle<std::vector<pat::MET> >& mets,
0039 std::vector<int>& match,
0040 const unsigned int iComb) {
0041
0042
0043
0044 for (unsigned idx = 0; idx < match.size(); ++idx) {
0045 if (isValid(match[idx], jets)) {
0046 switch (idx) {
0047 case TtFullLepEvtPartons::B:
0048 b_ = makeCandidate(jets, match[idx], jetCorrectionLevel_);
0049 break;
0050 case TtFullLepEvtPartons::BBar:
0051 bBar_ = makeCandidate(jets, match[idx], jetCorrectionLevel_);
0052 break;
0053 }
0054 }
0055 }
0056
0057
0058
0059
0060
0061 edm::Handle<TtGenEvent> genEvt;
0062 evt.getByToken(genEvtToken_, genEvt);
0063
0064
0065 if (!genEvt->isFullLeptonic() || !genEvt->lepton() || !genEvt->leptonBar()) {
0066 match.push_back(-1);
0067 match.push_back(-1);
0068 match.push_back(-1);
0069 match.push_back(-1);
0070 } else if (genEvt->isFullLeptonic(WDecay::kElec, WDecay::kElec) && elecs->size() >= 2) {
0071
0072 int iLepBar = findMatchingLepton(genEvt->leptonBar(), elecs);
0073 leptonBar_ = makeCandidate(elecs, iLepBar);
0074 match.push_back(iLepBar);
0075 int iLep = findMatchingLepton(genEvt->lepton(), elecs);
0076 lepton_ = makeCandidate(elecs, iLep);
0077 match.push_back(iLep);
0078
0079
0080 match.push_back(-1);
0081 match.push_back(-1);
0082 } else if (genEvt->isFullLeptonic(WDecay::kElec, WDecay::kMuon) && !elecs->empty() && !mus->empty()) {
0083 if (genEvt->leptonBar()->isElectron()) {
0084
0085 int iLepBar = findMatchingLepton(genEvt->leptonBar(), elecs);
0086 leptonBar_ = makeCandidate(elecs, iLepBar);
0087 match.push_back(iLepBar);
0088
0089 match.push_back(-1);
0090 match.push_back(-1);
0091
0092 int iLep = findMatchingLepton(genEvt->lepton(), mus);
0093 lepton_ = makeCandidate(mus, iLep);
0094 match.push_back(iLep);
0095 } else {
0096
0097 match.push_back(-1);
0098
0099 int iLepBar = findMatchingLepton(genEvt->leptonBar(), mus);
0100 leptonBar_ = makeCandidate(mus, iLepBar);
0101 match.push_back(iLepBar);
0102
0103 int iLep = findMatchingLepton(genEvt->lepton(), elecs);
0104 lepton_ = makeCandidate(elecs, iLep);
0105 match.push_back(iLep);
0106
0107 match.push_back(-1);
0108 }
0109 } else if (genEvt->isFullLeptonic(WDecay::kMuon, WDecay::kMuon) && mus->size() >= 2) {
0110
0111 match.push_back(-1);
0112 match.push_back(-1);
0113
0114
0115 int iLepBar = findMatchingLepton(genEvt->leptonBar(), mus);
0116 leptonBar_ = makeCandidate(mus, iLepBar);
0117 match.push_back(iLepBar);
0118 int iLep = findMatchingLepton(genEvt->lepton(), mus);
0119 lepton_ = makeCandidate(mus, iLep);
0120 match.push_back(iLep);
0121 } else {
0122 match.push_back(-1);
0123 match.push_back(-1);
0124 match.push_back(-1);
0125 match.push_back(-1);
0126 }
0127
0128
0129
0130
0131 if (!mets->empty()) {
0132
0133 buildMatchingNeutrinos(evt, mets);
0134 }
0135 }
0136
0137 template <typename O>
0138 int TtFullLepHypGenMatch::findMatchingLepton(const reco::GenParticle* genLep,
0139 const edm::Handle<std::vector<O> >& leps) {
0140 int idx = -1;
0141 double minDR = -1;
0142 for (unsigned i = 0; i < leps->size(); ++i) {
0143 double dR = deltaR(genLep->eta(), genLep->phi(), (*leps)[i].eta(), (*leps)[i].phi());
0144 if (minDR < 0 || dR < minDR) {
0145 minDR = dR;
0146 idx = i;
0147 }
0148 }
0149 return idx;
0150 }
0151
0152 void TtFullLepHypGenMatch::buildMatchingNeutrinos(edm::Event& evt, const edm::Handle<std::vector<pat::MET> >& mets) {
0153
0154 edm::Handle<TtGenEvent> genEvt;
0155 evt.getByToken(genEvtToken_, genEvt);
0156
0157 if (genEvt->isTtBar() && genEvt->isFullLeptonic() && genEvt->neutrino() && genEvt->neutrinoBar()) {
0158 double momXNu = genEvt->neutrino()->px();
0159 double momYNu = genEvt->neutrino()->py();
0160 double momXNuBar = genEvt->neutrinoBar()->px();
0161 double momYNuBar = genEvt->neutrinoBar()->py();
0162
0163 double momXMet = mets->at(0).px();
0164 double momYMet = mets->at(0).py();
0165
0166 double momXNeutrino = 0.5 * (momXNu - momXNuBar + momXMet);
0167 double momYNeutrino = 0.5 * (momYNu - momYNuBar + momYMet);
0168 double momXNeutrinoBar = momXMet - momXNeutrino;
0169 double momYNeutrinoBar = momYMet - momYNeutrino;
0170
0171 math::XYZTLorentzVector recNuFM(
0172 momXNeutrino, momYNeutrino, 0, sqrt(momXNeutrino * momXNeutrino + momYNeutrino * momYNeutrino));
0173 recNu = std::make_unique<reco::LeafCandidate>(0, recNuFM);
0174
0175 math::XYZTLorentzVector recNuBarFM(momXNeutrinoBar,
0176 momYNeutrinoBar,
0177 0,
0178 sqrt(momXNeutrinoBar * momXNeutrinoBar + momYNeutrinoBar * momYNeutrinoBar));
0179 recNuBar = std::make_unique<reco::LeafCandidate>(0, recNuBarFM);
0180 }
0181 }
0182
0183 #include "FWCore/Framework/interface/MakerMacros.h"
0184 DEFINE_FWK_MODULE(TtFullLepHypGenMatch);