Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef TtFullLepHypothesis_h
0002 #define TtFullLepHypothesis_h
0003 
0004 #include <memory>
0005 #include <vector>
0006 
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/stream/EDProducer.h"
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 
0013 #include "DataFormats/PatCandidates/interface/Electron.h"
0014 #include "DataFormats/PatCandidates/interface/Muon.h"
0015 #include "DataFormats/PatCandidates/interface/Jet.h"
0016 #include "DataFormats/PatCandidates/interface/MET.h"
0017 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0018 #include "DataFormats/Candidate/interface/ShallowClonePtrCandidate.h"
0019 
0020 #include "AnalysisDataFormats/TopObjects/interface/TtFullLeptonicEvent.h"
0021 
0022 /*
0023    \class   TtFullLepHypothesis TtFullLepHypothesis.h "TopQuarkAnalysis/TopJetCombination/interface/TtFullLepHypothesis.h"
0024 
0025    \brief   Interface class for the creation of full leptonic ttbar event hypotheses
0026 
0027    The class provides an interface for the creation of full leptonic ttbar event hypotheses. Input information is read
0028    from the event content and the proper candidate creation is taken care of. Hypotheses are characterized by the
0029    CompositeCandidate made of a ttbar pair (including all its decay products in a parton level interpretation) and an
0030    enumerator type key to specify the algorithm to determine the candidate (hypothesis class). The buildKey and the
0031    buildHypo class have to implemented by derived classes.
0032 **/
0033 
0034 class TtFullLepHypothesis : public edm::stream::EDProducer<> {
0035 public:
0036   /// default constructor
0037   explicit TtFullLepHypothesis(const edm::ParameterSet&);
0038 
0039 protected:
0040   /// produce the event hypothesis as CompositeCandidate and Key
0041   void produce(edm::Event&, const edm::EventSetup&) override;
0042   /// reset candidate pointers before hypo build process
0043   void resetCandidates();
0044   /// use one object in a collection to set a ShallowClonePtrCandidate
0045   template <typename C>
0046   std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<C>& handle, const int& idx);
0047   /// use one object in a jet collection to set a ShallowClonePtrCandidate with proper jet corrections
0048   std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<std::vector<pat::Jet> >& handle,
0049                                                                 const int& idx,
0050                                                                 const std::string& correctionLevel);
0051   /// return key
0052   int key() const { return key_; };
0053   /// return event hypothesis
0054   reco::CompositeCandidate hypo();
0055   /// check if index is in valid range of selected jets
0056   bool isValid(const int& idx, const edm::Handle<std::vector<pat::Jet> >& jets) {
0057     return (0 <= idx && idx < (int)jets->size());
0058   };
0059 
0060   // -----------------------------------------
0061   // implemet the following two functions
0062   // for a concrete event hypothesis
0063   // -----------------------------------------
0064 
0065   /// build the event hypothesis key
0066   virtual void buildKey() = 0;
0067   /// build event hypothesis from the reco objects of a semi-leptonic event
0068   virtual void buildHypo(edm::Event& evt,
0069                          const edm::Handle<std::vector<pat::Electron> >& elecs,
0070                          const edm::Handle<std::vector<pat::Muon> >& mus,
0071                          const edm::Handle<std::vector<pat::Jet> >& jets,
0072                          const edm::Handle<std::vector<pat::MET> >& mets,
0073                          std::vector<int>& match,
0074                          const unsigned int iComb) = 0;
0075 
0076 protected:
0077   /// internal check whether the match information exists or not,
0078   /// if false a blind dummy match vector will be used internally
0079   bool getMatch_;
0080   /// input label for all necessary collections
0081   edm::EDGetTokenT<std::vector<std::vector<int> > > matchToken_;
0082   edm::EDGetTokenT<std::vector<pat::Electron> > elecsToken_;
0083   edm::EDGetTokenT<std::vector<pat::Muon> > musToken_;
0084   edm::EDGetTokenT<std::vector<pat::Jet> > jetsToken_;
0085   edm::EDGetTokenT<std::vector<pat::MET> > metsToken_;
0086   /// specify the desired jet correction level (the default should
0087   /// be L3Absolute-'abs')
0088   std::string jetCorrectionLevel_;
0089   /// hypothesis key (to be set by the buildKey function)
0090   int key_;
0091   /// candidates for internal use for the creation of the hypothesis
0092   /// candidate
0093   std::unique_ptr<reco::ShallowClonePtrCandidate> lepton_;
0094   std::unique_ptr<reco::ShallowClonePtrCandidate> leptonBar_;
0095   std::unique_ptr<reco::ShallowClonePtrCandidate> b_;
0096   std::unique_ptr<reco::ShallowClonePtrCandidate> bBar_;
0097   std::unique_ptr<reco::ShallowClonePtrCandidate> neutrino_;
0098   std::unique_ptr<reco::ShallowClonePtrCandidate> neutrinoBar_;
0099   //reco::ShallowClonePtrCandidate *met_;
0100 
0101   /// candidates needed for the genmatch hypothesis
0102   std::unique_ptr<reco::LeafCandidate> recNu;
0103   std::unique_ptr<reco::LeafCandidate> recNuBar;
0104 };
0105 
0106 // unfortunately this has to be placed in the header since otherwise the function template
0107 // would cause unresolved references in classes derived from this base class
0108 template <typename C>
0109 std::unique_ptr<reco::ShallowClonePtrCandidate> TtFullLepHypothesis::makeCandidate(const edm::Handle<C>& handle,
0110                                                                                    const int& idx) {
0111   typedef typename C::value_type O;
0112   edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
0113   return std::make_unique<reco::ShallowClonePtrCandidate>(ptr, ptr->charge(), ptr->p4(), ptr->vertex());
0114 }
0115 
0116 #endif