Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef TtSemiLepHypothesis_h
0002 #define TtSemiLepHypothesis_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/PatCandidates/interface/MET.h"
0014 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0015 #include "DataFormats/Candidate/interface/ShallowClonePtrCandidate.h"
0016 
0017 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLeptonicEvent.h"
0018 
0019 /*
0020    \class   TtSemiLepHypothesis TtSemiLepHypothesis.h "TopQuarkAnalysis/TopJetCombination/interface/TtSemiLepHypothesis.h"
0021 
0022    \brief   Interface class for the creation of semi-leptonic ttbar event hypotheses
0023 
0024    The class provides an interface for the creation of semi-leptonic ttbar event hypotheses. Input information is read
0025    from the event content and the proper candidate creation is taken care of. Hypotheses are characterized by the
0026    CompositeCandidate made of a ttbar pair (including all its decay products in a parton level interpretation) and an
0027    enumerator type key to specify the algorithm that was used to determine the candidate (the "hypothesis class").
0028    The buildKey and the buildHypo methods have to implemented by derived classes.
0029 **/
0030 
0031 class TtSemiLepHypothesis : public edm::stream::EDProducer<> {
0032 public:
0033   /// default constructor
0034   explicit TtSemiLepHypothesis(const edm::ParameterSet&);
0035 
0036 protected:
0037   /// produce the event hypothesis as CompositeCandidate and Key
0038   void produce(edm::Event&, const edm::EventSetup&) override;
0039   /// reset candidate pointers before hypo build process
0040   void resetCandidates();
0041   /// helper function to construct the proper correction level string for corresponding quarkType,
0042   /// for unknown quarkTypes an empty string is returned
0043   std::string jetCorrectionLevel(const std::string& quarkType);
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   /// set neutrino, using mW = 80.4 to calculate the neutrino pz
0052   void setNeutrino(const edm::Handle<std::vector<pat::MET> >& met,
0053                    const edm::Handle<edm::View<reco::RecoCandidate> >& leps,
0054                    const int& idx,
0055                    const int& type);
0056   /// return key
0057   /// minimalistic build function for simple hypotheses
0058   void buildHypo(const edm::Handle<edm::View<reco::RecoCandidate> >& leps,
0059                  const edm::Handle<std::vector<pat::MET> >& mets,
0060                  const edm::Handle<std::vector<pat::Jet> >& jets,
0061                  std::vector<int>& jetPartonAssociation);
0062   int key() const { return key_; };
0063   /// return event hypothesis
0064   reco::CompositeCandidate hypo();
0065   /// check if index is in valid range of selected jets
0066   bool isValid(const int& idx, const edm::Handle<std::vector<pat::Jet> >& jets) {
0067     return (0 <= idx && idx < (int)jets->size());
0068   };
0069   /// determine lepton type of reco candidate and return a corresponding WDecay::LeptonType; the type is kNone if it is whether a muon nor an electron
0070   WDecay::LeptonType leptonType(const reco::RecoCandidate* cand);
0071 
0072   // -----------------------------------------
0073   // implemet the following two functions
0074   // for a concrete event hypothesis
0075   // -----------------------------------------
0076 
0077   /// build the event hypothesis key
0078   virtual void buildKey() = 0;
0079   /// build event hypothesis from the reco objects of a semi-leptonic event
0080   virtual void buildHypo(edm::Event& event,
0081                          const edm::Handle<edm::View<reco::RecoCandidate> >& lepton,
0082                          const edm::Handle<std::vector<pat::MET> >& neutrino,
0083                          const edm::Handle<std::vector<pat::Jet> >& jets,
0084                          std::vector<int>& jetPartonAssociation,
0085                          const unsigned int iComb) = 0;
0086 
0087 protected:
0088   /// internal check whether the match information exists or not,
0089   /// if false a blind dummy match vector will be used internally
0090   bool getMatch_;
0091   /// input label for all necessary collections
0092   edm::EDGetTokenT<std::vector<pat::Jet> > jetsToken_;
0093   edm::EDGetTokenT<edm::View<reco::RecoCandidate> > lepsToken_;
0094   edm::EDGetTokenT<std::vector<pat::MET> > metsToken_;
0095   edm::EDGetTokenT<std::vector<std::vector<int> > > matchToken_;
0096   edm::EDGetTokenT<int> nJetsConsideredToken_;
0097   /// specify the desired jet correction level (the default should
0098   /// be L3Absolute-'abs')
0099   std::string jetCorrectionLevel_;
0100   /// hypothesis key (to be set by the buildKey function)
0101   int key_;
0102   /// algorithm used to calculate neutrino solutions (see cfi for further details)
0103   int neutrinoSolutionType_;
0104   /// number of real neutrino solutions:
0105   /// -1 if not determined, 0 if only complex, 2 if two real solutions
0106   int numberOfRealNeutrinoSolutions_;
0107   /// candidates for internal use for the creation of the hypothesis
0108   /// candidate
0109   std::unique_ptr<reco::ShallowClonePtrCandidate> lightQ_;
0110   std::unique_ptr<reco::ShallowClonePtrCandidate> lightQBar_;
0111   std::unique_ptr<reco::ShallowClonePtrCandidate> hadronicB_;
0112   std::unique_ptr<reco::ShallowClonePtrCandidate> leptonicB_;
0113   std::unique_ptr<reco::ShallowClonePtrCandidate> neutrino_;
0114   std::unique_ptr<reco::ShallowClonePtrCandidate> lepton_;
0115 };
0116 
0117 // has to be placed in the header since otherwise the function template
0118 // would cause unresolved references in classes derived from this base class
0119 template <typename C>
0120 std::unique_ptr<reco::ShallowClonePtrCandidate> TtSemiLepHypothesis::makeCandidate(const edm::Handle<C>& handle,
0121                                                                                    const int& idx) {
0122   typedef typename C::value_type O;
0123   edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
0124   return std::make_unique<reco::ShallowClonePtrCandidate>(ptr, ptr->charge(), ptr->p4(), ptr->vertex());
0125 }
0126 
0127 #endif