CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtSemiLepHypothesis.h
Go to the documentation of this file.
1 #ifndef TtSemiLepHypothesis_h
2 #define TtSemiLepHypothesis_h
3 
4 #include <memory>
5 #include <vector>
6 
11 
16 
18 
19 /*
20  \class TtSemiLepHypothesis TtSemiLepHypothesis.h "TopQuarkAnalysis/TopJetCombination/interface/TtSemiLepHypothesis.h"
21 
22  \brief Interface class for the creation of semi-leptonic ttbar event hypotheses
23 
24  The class provides an interface for the creation of semi-leptonic ttbar event hypotheses. Input information is read
25  from the event content and the proper candidate creation is taken care of. Hypotheses are characterized by the
26  CompositeCandidate made of a ttbar pair (including all its decay products in a parton level interpretation) and an
27  enumerator type key to specify the algorithm that was used to determine the candidate (the "hypothesis class").
28  The buildKey and the buildHypo methods have to implemented by derived classes.
29 **/
30 
32 public:
34  explicit TtSemiLepHypothesis(const edm::ParameterSet&);
35 
36 protected:
38  void produce(edm::Event&, const edm::EventSetup&) override;
40  void resetCandidates();
43  std::string jetCorrectionLevel(const std::string& quarkType);
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  void setNeutrino(const edm::Handle<std::vector<pat::MET> >& met,
54  const int& idx,
55  const int& type);
59  const edm::Handle<std::vector<pat::MET> >& mets,
60  const edm::Handle<std::vector<pat::Jet> >& jets,
61  std::vector<int>& jetPartonAssociation);
62  int key() const { return key_; };
66  bool isValid(const int& idx, const edm::Handle<std::vector<pat::Jet> >& jets) {
67  return (0 <= idx && idx < (int)jets->size());
68  };
71 
72  // -----------------------------------------
73  // implemet the following two functions
74  // for a concrete event hypothesis
75  // -----------------------------------------
76 
78  virtual void buildKey() = 0;
80  virtual void buildHypo(edm::Event& event,
82  const edm::Handle<std::vector<pat::MET> >& neutrino,
83  const edm::Handle<std::vector<pat::Jet> >& jets,
84  std::vector<int>& jetPartonAssociation,
85  const unsigned int iComb) = 0;
86 
87 protected:
90  bool getMatch_;
101  int key_;
109  std::unique_ptr<reco::ShallowClonePtrCandidate> lightQ_;
110  std::unique_ptr<reco::ShallowClonePtrCandidate> lightQBar_;
111  std::unique_ptr<reco::ShallowClonePtrCandidate> hadronicB_;
112  std::unique_ptr<reco::ShallowClonePtrCandidate> leptonicB_;
113  std::unique_ptr<reco::ShallowClonePtrCandidate> neutrino_;
114  std::unique_ptr<reco::ShallowClonePtrCandidate> lepton_;
115 };
116 
117 // has to be placed in the header since otherwise the function template
118 // would cause unresolved references in classes derived from this base class
119 template <typename C>
120 std::unique_ptr<reco::ShallowClonePtrCandidate> TtSemiLepHypothesis::makeCandidate(const edm::Handle<C>& handle,
121  const int& idx) {
122  typedef typename C::value_type O;
123  edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
124  return std::make_unique<reco::ShallowClonePtrCandidate>(ptr, ptr->charge(), ptr->p4(), ptr->vertex());
125 }
126 
127 #endif
std::string jetCorrectionLevel_
edm::EDGetTokenT< int > nJetsConsideredToken_
WDecay::LeptonType leptonType(const reco::RecoCandidate *cand)
determine lepton type of reco candidate and return a corresponding WDecay::LeptonType; the type is kN...
bool isValid(const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
check if index is in valid range of selected jets
void buildHypo(const edm::Handle< edm::View< reco::RecoCandidate > > &leps, const edm::Handle< std::vector< pat::MET > > &mets, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &jetPartonAssociation)
minimalistic build function for simple hypotheses
std::unique_ptr< reco::ShallowClonePtrCandidate > leptonicB_
TtSemiLepHypothesis(const edm::ParameterSet &)
default constructor
std::unique_ptr< reco::ShallowClonePtrCandidate > lepton_
void resetCandidates()
reset candidate pointers before hypo build process
vector< PseudoJet > jets
std::unique_ptr< reco::ShallowClonePtrCandidate > makeCandidate(const edm::Handle< C > &handle, const int &idx)
use one object in a collection to set a ShallowClonePtrCandidate
tuple handle
Definition: patZpeak.py:25
virtual void buildKey()=0
build the event hypothesis key
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
input label for all necessary collections
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQBar_
void setNeutrino(const edm::Handle< std::vector< pat::MET > > &met, const edm::Handle< edm::View< reco::RecoCandidate > > &leps, const int &idx, const int &type)
set neutrino, using mW = 80.4 to calculate the neutrino pz
reco::CompositeCandidate hypo()
return event hypothesis
std::unique_ptr< reco::ShallowClonePtrCandidate > neutrino_
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQ_
edm::EDGetTokenT< edm::View< reco::RecoCandidate > > lepsToken_
edm::EDGetTokenT< std::vector< std::vector< int > > > matchToken_
int neutrinoSolutionType_
algorithm used to calculate neutrino solutions (see cfi for further details)
void produce(edm::Event &, const edm::EventSetup &) override
produce the event hypothesis as CompositeCandidate and Key
std::unique_ptr< reco::ShallowClonePtrCandidate > hadronicB_
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...
int key_
hypothesis key (to be set by the buildKey function)