Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "TopQuarkAnalysis/TopJetCombination/interface/TtFullHadHypothesis.h"
0002 
0003 #include "DataFormats/PatCandidates/interface/Particle.h"
0004 
0005 class TtFullHadHypKinFit : public TtFullHadHypothesis {
0006 public:
0007   explicit TtFullHadHypKinFit(const edm::ParameterSet&);
0008 
0009 private:
0010   /// build the event hypothesis key
0011   void buildKey() override { key_ = TtFullHadronicEvent::kKinFit; };
0012   /// build event hypothesis from the reco objects of a full-hadronic event
0013   void buildHypo(edm::Event&,
0014                  const edm::Handle<std::vector<pat::Jet> >&,
0015                  std::vector<int>&,
0016                  const unsigned int iComb) override;
0017 
0018   edm::EDGetTokenT<std::vector<int> > statusToken_;
0019   edm::EDGetTokenT<std::vector<pat::Particle> > lightQToken_;
0020   edm::EDGetTokenT<std::vector<pat::Particle> > lightQBarToken_;
0021   edm::EDGetTokenT<std::vector<pat::Particle> > bToken_;
0022   edm::EDGetTokenT<std::vector<pat::Particle> > bBarToken_;
0023   edm::EDGetTokenT<std::vector<pat::Particle> > lightPToken_;
0024   edm::EDGetTokenT<std::vector<pat::Particle> > lightPBarToken_;
0025 };
0026 
0027 TtFullHadHypKinFit::TtFullHadHypKinFit(const edm::ParameterSet& cfg)
0028     : TtFullHadHypothesis(cfg),
0029       statusToken_(consumes<std::vector<int> >(cfg.getParameter<edm::InputTag>("status"))),
0030       lightQToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightQTag"))),
0031       lightQBarToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightQBarTag"))),
0032       bToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("bTag"))),
0033       bBarToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("bBarTag"))),
0034       lightPToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightPTag"))),
0035       lightPBarToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightPBarTag"))) {}
0036 
0037 void TtFullHadHypKinFit::buildHypo(edm::Event& evt,
0038                                    const edm::Handle<std::vector<pat::Jet> >& jets,
0039                                    std::vector<int>& match,
0040                                    const unsigned int iComb) {
0041   edm::Handle<std::vector<int> > status;
0042   evt.getByToken(statusToken_, status);
0043   if ((*status)[iComb] != 0) {
0044     // create empty hypothesis if kinematic fit did not converge
0045     return;
0046   }
0047 
0048   edm::Handle<std::vector<pat::Particle> > lightQ;
0049   edm::Handle<std::vector<pat::Particle> > lightQBar;
0050   edm::Handle<std::vector<pat::Particle> > b;
0051   edm::Handle<std::vector<pat::Particle> > bBar;
0052   edm::Handle<std::vector<pat::Particle> > lightP;
0053   edm::Handle<std::vector<pat::Particle> > lightPBar;
0054 
0055   evt.getByToken(lightQToken_, lightQ);
0056   evt.getByToken(lightQBarToken_, lightQBar);
0057   evt.getByToken(bToken_, b);
0058   evt.getByToken(bBarToken_, bBar);
0059   evt.getByToken(lightPToken_, lightP);
0060   evt.getByToken(lightPBarToken_, lightPBar);
0061 
0062   // -----------------------------------------------------
0063   // add jets
0064   // -----------------------------------------------------
0065   if (!(lightQ->empty() || lightQBar->empty() || b->empty() || bBar->empty() || lightP->empty() ||
0066         lightPBar->empty())) {
0067     lightQ_ = makeCandidate(lightQ, iComb);
0068     lightQBar_ = makeCandidate(lightQBar, iComb);
0069     b_ = makeCandidate(b, iComb);
0070     bBar_ = makeCandidate(bBar, iComb);
0071     lightP_ = makeCandidate(lightP, iComb);
0072     lightPBar_ = makeCandidate(lightPBar, iComb);
0073   }
0074 }
0075 
0076 #include "FWCore/Framework/interface/MakerMacros.h"
0077 DEFINE_FWK_MODULE(TtFullHadHypKinFit);