Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:20

0001 #include "DataFormats/Math/interface/deltaR.h"
0002 #include "AnalysisDataFormats/TopObjects/interface/TtGenEvent.h"
0003 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
0004 #include "TopQuarkAnalysis/TopJetCombination/interface/TtSemiLepHypothesis.h"
0005 
0006 class TtSemiLepHypGenMatch : public TtSemiLepHypothesis {
0007 public:
0008   explicit TtSemiLepHypGenMatch(const edm::ParameterSet&);
0009 
0010 private:
0011   /// build the event hypothesis key
0012   void buildKey() override { key_ = TtSemiLeptonicEvent::kGenMatch; };
0013   /// build event hypothesis from the reco objects of a semi-leptonic event
0014   void buildHypo(edm::Event&,
0015                  const edm::Handle<edm::View<reco::RecoCandidate> >&,
0016                  const edm::Handle<std::vector<pat::MET> >&,
0017                  const edm::Handle<std::vector<pat::Jet> >&,
0018                  std::vector<int>&,
0019                  const unsigned int iComb) override;
0020   /// find index of the candidate nearest to the singleLepton of the generator event in the collection; return -1 if this fails
0021   int findMatchingLepton(const edm::Handle<TtGenEvent>& genEvt, const edm::Handle<edm::View<reco::RecoCandidate> >&);
0022 
0023 protected:
0024   edm::EDGetTokenT<TtGenEvent> genEvtToken_;
0025 };
0026 
0027 TtSemiLepHypGenMatch::TtSemiLepHypGenMatch(const edm::ParameterSet& cfg)
0028     : TtSemiLepHypothesis(cfg), genEvtToken_(consumes<TtGenEvent>(edm::InputTag("genEvt"))) {}
0029 
0030 void TtSemiLepHypGenMatch::buildHypo(edm::Event& evt,
0031                                      const edm::Handle<edm::View<reco::RecoCandidate> >& leps,
0032                                      const edm::Handle<std::vector<pat::MET> >& mets,
0033                                      const edm::Handle<std::vector<pat::Jet> >& jets,
0034                                      std::vector<int>& match,
0035                                      const unsigned int iComb) {
0036   // -----------------------------------------------------
0037   // get genEvent (to distinguish between uds and c quarks
0038   // and for the lepton matching)
0039   // -----------------------------------------------------
0040   edm::Handle<TtGenEvent> genEvt;
0041   evt.getByToken(genEvtToken_, genEvt);
0042 
0043   // -----------------------------------------------------
0044   // add jets
0045   // -----------------------------------------------------
0046   for (unsigned idx = 0; idx < match.size(); ++idx) {
0047     if (isValid(match[idx], jets)) {
0048       switch (idx) {
0049         case TtSemiLepEvtPartons::LightQ:
0050           if (std::abs(genEvt->hadronicDecayQuark()->pdgId()) == 4)
0051             lightQ_ = makeCandidate(jets, match[idx], jetCorrectionLevel("cQuark"));
0052           else
0053             lightQ_ = makeCandidate(jets, match[idx], jetCorrectionLevel("udsQuark"));
0054           break;
0055         case TtSemiLepEvtPartons::LightQBar:
0056           if (std::abs(genEvt->hadronicDecayQuarkBar()->pdgId()) == 4)
0057             lightQBar_ = makeCandidate(jets, match[idx], jetCorrectionLevel("cQuark"));
0058           else
0059             lightQBar_ = makeCandidate(jets, match[idx], jetCorrectionLevel("udsQuark"));
0060           break;
0061         case TtSemiLepEvtPartons::HadB:
0062           hadronicB_ = makeCandidate(jets, match[idx], jetCorrectionLevel("bQuark"));
0063           break;
0064         case TtSemiLepEvtPartons::LepB:
0065           leptonicB_ = makeCandidate(jets, match[idx], jetCorrectionLevel("bQuark"));
0066           break;
0067       }
0068     }
0069   }
0070 
0071   // -----------------------------------------------------
0072   // add lepton
0073   // -----------------------------------------------------
0074   int iLepton = findMatchingLepton(genEvt, leps);
0075   if (iLepton < 0)
0076     return;
0077   lepton_ = makeCandidate(leps, iLepton);
0078   match.push_back(iLepton);
0079 
0080   // -----------------------------------------------------
0081   // add neutrino
0082   // -----------------------------------------------------
0083   if (mets->empty())
0084     return;
0085   if (neutrinoSolutionType_ == -1)
0086     neutrino_ = makeCandidate(mets, 0);
0087   else
0088     setNeutrino(mets, leps, iLepton, neutrinoSolutionType_);
0089 }
0090 
0091 /// find index of the candidate nearest to the singleLepton of the generator event in the collection; return -1 if this fails
0092 int TtSemiLepHypGenMatch::findMatchingLepton(const edm::Handle<TtGenEvent>& genEvt,
0093                                              const edm::Handle<edm::View<reco::RecoCandidate> >& leps) {
0094   int genIdx = -1;
0095 
0096   // jump out with -1 when the collection is empty
0097   if (leps->empty())
0098     return genIdx;
0099 
0100   if (genEvt->isTtBar() && genEvt->isSemiLeptonic(leptonType(&(leps->front()))) && genEvt->singleLepton()) {
0101     double minDR = -1;
0102     for (unsigned i = 0; i < leps->size(); ++i) {
0103       double dR =
0104           deltaR(genEvt->singleLepton()->eta(), genEvt->singleLepton()->phi(), (*leps)[i].eta(), (*leps)[i].phi());
0105       if (minDR < 0 || dR < minDR) {
0106         minDR = dR;
0107         genIdx = i;
0108       }
0109     }
0110   }
0111   return genIdx;
0112 }
0113 
0114 #include "FWCore/Framework/interface/MakerMacros.h"
0115 DEFINE_FWK_MODULE(TtSemiLepHypGenMatch);