Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef TtFullHadHypothesis_h
0002 #define TtFullHadHypothesis_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 
0012 #include "DataFormats/PatCandidates/interface/Jet.h"
0013 #include "DataFormats/Candidate/interface/ShallowClonePtrCandidate.h"
0014 
0015 #include "AnalysisDataFormats/TopObjects/interface/TtFullHadronicEvent.h"
0016 
0017 /*
0018    \class   TtFullHadHypothesis TtFullHadHypothesis.h "TopQuarkAnalysis/TopJetCombination/interface/TtFullHadHypothesis.h"
0019 
0020    \brief   Interface class for the creation of full-hadronic ttbar event hypotheses
0021 
0022    The class provides an interface for the creation of full-hadronic ttbar event hypotheses. Input information is read
0023    from the event content and the proper candidate creation is taken care of. Hypotheses are characterized by the
0024    CompositeCandidate made of a ttbar pair (including all its decay products in a parton level interpretation) and an
0025    enumerator type key to specify the algorithm to determine the candidate (hypothesis cklass). The buildKey and the
0026    buildHypo class have to implemented by derived classes.
0027 **/
0028 
0029 class TtFullHadHypothesis : public edm::stream::EDProducer<> {
0030 public:
0031   /// default constructor
0032   explicit TtFullHadHypothesis(const edm::ParameterSet& cfg);
0033 
0034 protected:
0035   /// produce the event hypothesis as CompositeCandidate and Key
0036   void produce(edm::Event&, const edm::EventSetup&) override;
0037   /// reset candidate pointers before hypo build process
0038   void resetCandidates();
0039   /// helper function to construct the proper correction level string for corresponding quarkType,
0040   /// for unknown quarkTypes an empty string is returned
0041   std::string jetCorrectionLevel(const std::string& quarkType);
0042   /// use one object in a collection to set a ShallowClonePtrCandidate
0043   template <typename C>
0044   std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<C>& handle, const int& idx);
0045   /// use one object in a jet collection to set a ShallowClonePtrCandidate with proper jet corrections
0046   std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<std::vector<pat::Jet> >& handle,
0047                                                                 const int& idx,
0048                                                                 const std::string& correctionLevel);
0049   /// return key
0050   int key() const { return key_; };
0051   /// return event hypothesis
0052   reco::CompositeCandidate hypo();
0053   /// check if index is in valid range of selected jets
0054   bool isValid(const int& idx, const edm::Handle<std::vector<pat::Jet> >& jets) {
0055     return (0 <= idx && idx < (int)jets->size());
0056   };
0057 
0058   // -----------------------------------------
0059   // implemet the following two functions
0060   // for a concrete event hypothesis
0061   // -----------------------------------------
0062 
0063   /// build the event hypothesis key
0064   virtual void buildKey() = 0;
0065   /// build event hypothesis from the reco objects of a full-hadronic event
0066   virtual void buildHypo(edm::Event& event,
0067                          const edm::Handle<std::vector<pat::Jet> >& jets,
0068                          std::vector<int>& jetPartonAssociation,
0069                          const unsigned int iComb) = 0;
0070 
0071 protected:
0072   /// internal check whether the match information exists or not,
0073   /// if false a blind dummy match vector will be used internally
0074   bool getMatch_;
0075   /// input label for all necessary collections
0076   edm::EDGetTokenT<std::vector<pat::Jet> > jetsToken_;
0077   edm::EDGetTokenT<std::vector<std::vector<int> > > matchToken_;
0078   /// specify the desired jet correction level (the default should
0079   /// be L3Absolute-'abs')
0080   std::string jetCorrectionLevel_;
0081   /// hypothesis key (to be set by the buildKey function)
0082   int key_;
0083   /// candidates for internal use for the creation of the hypothesis
0084   /// candidate
0085   std::unique_ptr<reco::ShallowClonePtrCandidate> lightQ_;
0086   std::unique_ptr<reco::ShallowClonePtrCandidate> lightQBar_;
0087   std::unique_ptr<reco::ShallowClonePtrCandidate> b_;
0088   std::unique_ptr<reco::ShallowClonePtrCandidate> bBar_;
0089   std::unique_ptr<reco::ShallowClonePtrCandidate> lightP_;
0090   std::unique_ptr<reco::ShallowClonePtrCandidate> lightPBar_;
0091 };
0092 
0093 // has to be placed in the header since otherwise the function template
0094 // would cause unresolved references in classes derived from this base class
0095 template <typename C>
0096 std::unique_ptr<reco::ShallowClonePtrCandidate> TtFullHadHypothesis::makeCandidate(const edm::Handle<C>& handle,
0097                                                                                    const int& idx) {
0098   typedef typename C::value_type O;
0099   edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
0100   return std::make_unique<reco::ShallowClonePtrCandidate>(ptr, ptr->charge(), ptr->p4(), ptr->vertex());
0101 }
0102 #endif