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&);
36  ~TtSemiLepHypothesis() override;
37 
38 protected:
40  void produce(edm::Event&, const edm::EventSetup&) override;
42  void resetCandidates();
45  std::string jetCorrectionLevel(const std::string& quarkType);
47  template <typename C>
50  void setCandidate(const edm::Handle<std::vector<pat::Jet> >& handle,
51  const int& idx,
53  const std::string& correctionLevel);
55  void setNeutrino(const edm::Handle<std::vector<pat::MET> >& met,
57  const int& idx,
58  const int& type);
62  const edm::Handle<std::vector<pat::MET> >& mets,
63  const edm::Handle<std::vector<pat::Jet> >& jets,
64  std::vector<int>& jetPartonAssociation);
65  int key() const { return key_; };
69  bool isValid(const int& idx, const edm::Handle<std::vector<pat::Jet> >& jets) {
70  return (0 <= idx && idx < (int)jets->size());
71  };
74 
75  // -----------------------------------------
76  // implemet the following two functions
77  // for a concrete event hypothesis
78  // -----------------------------------------
79 
81  virtual void buildKey() = 0;
83  virtual void buildHypo(edm::Event& event,
85  const edm::Handle<std::vector<pat::MET> >& neutrino,
86  const edm::Handle<std::vector<pat::Jet> >& jets,
87  std::vector<int>& jetPartonAssociation,
88  const unsigned int iComb) = 0;
89 
90 protected:
93  bool getMatch_;
104  int key_;
118 };
119 
120 // has to be placed in the header since otherwise the function template
121 // would cause unresolved references in classes derived from this base class
122 template <typename C>
124  const int& idx,
126  typedef typename C::value_type O;
127  edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
128  clone = new reco::ShallowClonePtrCandidate(ptr, ptr->charge(), ptr->p4(), ptr->vertex());
129 }
130 
131 #endif
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
reco::ShallowClonePtrCandidate * lepton_
reco::ShallowClonePtrCandidate * lightQBar_
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
TtSemiLepHypothesis(const edm::ParameterSet &)
default constructor
reco::ShallowClonePtrCandidate * neutrino_
void resetCandidates()
reset candidate pointers before hypo build process
vector< PseudoJet > jets
tuple handle
Definition: patZpeak.py:23
virtual void buildKey()=0
build the event hypothesis key
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
input label for all necessary collections
reco::ShallowClonePtrCandidate * hadronicB_
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
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
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
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
reco::ShallowClonePtrCandidate * lightQ_
~TtSemiLepHypothesis() override
default destructor
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
reco::ShallowClonePtrCandidate * leptonicB_
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)