CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtFullLepHypothesis.h
Go to the documentation of this file.
1 #ifndef TtFullLepHypothesis_h
2 #define TtFullLepHypothesis_h
3 
4 #include <memory>
5 #include <vector>
6 
12 
19 
21 
22 /*
23  \class TtFullLepHypothesis TtFullLepHypothesis.h "TopQuarkAnalysis/TopJetCombination/interface/TtFullLepHypothesis.h"
24 
25  \brief Interface class for the creation of full leptonic ttbar event hypotheses
26 
27  The class provides an interface for the creation of full leptonic ttbar event hypotheses. Input information is read
28  from the event content and the proper candidate creation is taken care of. Hypotheses are characterized by the
29  CompositeCandidate made of a ttbar pair (including all its decay products in a parton level interpretation) and an
30  enumerator type key to specify the algorithm to determine the candidate (hypothesis class). The buildKey and the
31  buildHypo class have to implemented by derived classes.
32 **/
33 
35 public:
37  explicit TtFullLepHypothesis(const edm::ParameterSet&);
38 
39 protected:
41  void produce(edm::Event&, const edm::EventSetup&) override;
43  void resetCandidates();
45  template <typename C>
46  std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<C>& handle, const int& idx);
48  std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<std::vector<pat::Jet> >& handle,
49  const int& idx,
50  const std::string& correctionLevel);
52  int key() const { return key_; };
56  bool isValid(const int& idx, const edm::Handle<std::vector<pat::Jet> >& jets) {
57  return (0 <= idx && idx < (int)jets->size());
58  };
59 
60  // -----------------------------------------
61  // implemet the following two functions
62  // for a concrete event hypothesis
63  // -----------------------------------------
64 
66  virtual void buildKey() = 0;
68  virtual void buildHypo(edm::Event& evt,
69  const edm::Handle<std::vector<pat::Electron> >& elecs,
70  const edm::Handle<std::vector<pat::Muon> >& mus,
71  const edm::Handle<std::vector<pat::Jet> >& jets,
72  const edm::Handle<std::vector<pat::MET> >& mets,
73  std::vector<int>& match,
74  const unsigned int iComb) = 0;
75 
76 protected:
79  bool getMatch_;
90  int key_;
93  std::unique_ptr<reco::ShallowClonePtrCandidate> lepton_;
94  std::unique_ptr<reco::ShallowClonePtrCandidate> leptonBar_;
95  std::unique_ptr<reco::ShallowClonePtrCandidate> b_;
96  std::unique_ptr<reco::ShallowClonePtrCandidate> bBar_;
97  std::unique_ptr<reco::ShallowClonePtrCandidate> neutrino_;
98  std::unique_ptr<reco::ShallowClonePtrCandidate> neutrinoBar_;
99  //reco::ShallowClonePtrCandidate *met_;
100 
102  std::unique_ptr<reco::LeafCandidate> recNu;
103  std::unique_ptr<reco::LeafCandidate> recNuBar;
104 };
105 
106 // unfortunately this has to be placed in the header since otherwise the function template
107 // would cause unresolved references in classes derived from this base class
108 template <typename C>
109 std::unique_ptr<reco::ShallowClonePtrCandidate> TtFullLepHypothesis::makeCandidate(const edm::Handle<C>& handle,
110  const int& idx) {
111  typedef typename C::value_type O;
112  edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
113  return std::make_unique<reco::ShallowClonePtrCandidate>(ptr, ptr->charge(), ptr->p4(), ptr->vertex());
114 }
115 
116 #endif
std::unique_ptr< reco::LeafCandidate > recNu
candidates needed for the genmatch hypothesis
std::unique_ptr< reco::ShallowClonePtrCandidate > leptonBar_
int key_
hypothesis key (to be set by the buildKey function)
virtual void buildKey()=0
build the event hypothesis key
vector< PseudoJet > jets
tuple handle
Definition: patZpeak.py:25
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
void resetCandidates()
reset candidate pointers before hypo build process
reco::CompositeCandidate hypo()
return event hypothesis
std::unique_ptr< reco::ShallowClonePtrCandidate > neutrino_
std::unique_ptr< reco::ShallowClonePtrCandidate > neutrinoBar_
edm::EDGetTokenT< std::vector< std::vector< int > > > matchToken_
input label for all necessary collections
virtual void buildHypo(edm::Event &evt, const edm::Handle< std::vector< pat::Electron > > &elecs, const edm::Handle< std::vector< pat::Muon > > &mus, const edm::Handle< std::vector< pat::Jet > > &jets, const edm::Handle< std::vector< pat::MET > > &mets, std::vector< int > &match, const unsigned int iComb)=0
build event hypothesis from the reco objects of a semi-leptonic event
edm::EDGetTokenT< std::vector< pat::Electron > > elecsToken_
void produce(edm::Event &, const edm::EventSetup &) override
produce the event hypothesis as CompositeCandidate and Key
std::unique_ptr< reco::LeafCandidate > recNuBar
std::string jetCorrectionLevel_
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
std::unique_ptr< reco::ShallowClonePtrCandidate > lepton_
int key() const
return key
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::unique_ptr< reco::ShallowClonePtrCandidate > b_
std::unique_ptr< reco::ShallowClonePtrCandidate > makeCandidate(const edm::Handle< C > &handle, const int &idx)
use one object in a collection to set a ShallowClonePtrCandidate
TtFullLepHypothesis(const edm::ParameterSet &)
default constructor
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets
edm::EDGetTokenT< std::vector< pat::Muon > > musToken_
std::unique_ptr< reco::ShallowClonePtrCandidate > bBar_