CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtFullHadHypothesis.h
Go to the documentation of this file.
1 #ifndef TtFullHadHypothesis_h
2 #define TtFullHadHypothesis_h
3 
4 #include <memory>
5 #include <vector>
6 
11 
14 
16 
17 /*
18  \class TtFullHadHypothesis TtFullHadHypothesis.h "TopQuarkAnalysis/TopJetCombination/interface/TtFullHadHypothesis.h"
19 
20  \brief Interface class for the creation of full-hadronic ttbar event hypotheses
21 
22  The class provides an interface for the creation of full-hadronic ttbar event hypotheses. Input information is read
23  from the event content and the proper candidate creation is taken care of. Hypotheses are characterized by the
24  CompositeCandidate made of a ttbar pair (including all its decay products in a parton level interpretation) and an
25  enumerator type key to specify the algorithm to determine the candidate (hypothesis cklass). The buildKey and the
26  buildHypo class have to implemented by derived classes.
27 **/
28 
30 public:
32  explicit TtFullHadHypothesis(const edm::ParameterSet& cfg);
33 
34 protected:
36  void produce(edm::Event&, const edm::EventSetup&) override;
38  void resetCandidates();
41  std::string jetCorrectionLevel(const std::string& quarkType);
43  template <typename C>
44  std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<C>& handle, const int& idx);
46  std::unique_ptr<reco::ShallowClonePtrCandidate> makeCandidate(const edm::Handle<std::vector<pat::Jet> >& handle,
47  const int& idx,
48  const std::string& correctionLevel);
50  int key() const { return key_; };
54  bool isValid(const int& idx, const edm::Handle<std::vector<pat::Jet> >& jets) {
55  return (0 <= idx && idx < (int)jets->size());
56  };
57 
58  // -----------------------------------------
59  // implemet the following two functions
60  // for a concrete event hypothesis
61  // -----------------------------------------
62 
64  virtual void buildKey() = 0;
66  virtual void buildHypo(edm::Event& event,
67  const edm::Handle<std::vector<pat::Jet> >& jets,
68  std::vector<int>& jetPartonAssociation,
69  const unsigned int iComb) = 0;
70 
71 protected:
74  bool getMatch_;
82  int key_;
85  std::unique_ptr<reco::ShallowClonePtrCandidate> lightQ_;
86  std::unique_ptr<reco::ShallowClonePtrCandidate> lightQBar_;
87  std::unique_ptr<reco::ShallowClonePtrCandidate> b_;
88  std::unique_ptr<reco::ShallowClonePtrCandidate> bBar_;
89  std::unique_ptr<reco::ShallowClonePtrCandidate> lightP_;
90  std::unique_ptr<reco::ShallowClonePtrCandidate> lightPBar_;
91 };
92 
93 // has to be placed in the header since otherwise the function template
94 // would cause unresolved references in classes derived from this base class
95 template <typename C>
96 std::unique_ptr<reco::ShallowClonePtrCandidate> TtFullHadHypothesis::makeCandidate(const edm::Handle<C>& handle,
97  const int& idx) {
98  typedef typename C::value_type O;
99  edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
100  return std::make_unique<reco::ShallowClonePtrCandidate>(ptr, ptr->charge(), ptr->p4(), ptr->vertex());
101 }
102 #endif
tuple cfg
Definition: looper.py:296
int key() const
return key
std::unique_ptr< reco::ShallowClonePtrCandidate > makeCandidate(const edm::Handle< C > &handle, const int &idx)
use one object in a collection to set a ShallowClonePtrCandidate
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...
std::unique_ptr< reco::ShallowClonePtrCandidate > b_
vector< PseudoJet > jets
std::unique_ptr< reco::ShallowClonePtrCandidate > lightP_
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets
tuple handle
Definition: patZpeak.py:25
virtual void buildHypo(edm::Event &event, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &jetPartonAssociation, const unsigned int iComb)=0
build event hypothesis from the reco objects of a full-hadronic event
std::unique_ptr< reco::ShallowClonePtrCandidate > bBar_
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQ_
void produce(edm::Event &, const edm::EventSetup &) override
produce the event hypothesis as CompositeCandidate and Key
edm::EDGetTokenT< std::vector< std::vector< int > > > matchToken_
std::string jetCorrectionLevel_
std::unique_ptr< reco::ShallowClonePtrCandidate > lightPBar_
virtual void buildKey()=0
build the event hypothesis key
void resetCandidates()
reset candidate pointers before hypo build process
int key_
hypothesis key (to be set by the buildKey function)
TtFullHadHypothesis(const edm::ParameterSet &cfg)
default constructor
reco::CompositeCandidate hypo()
return event hypothesis
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
input label for all necessary collections
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQBar_