CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtFullLepHypothesis.cc
Go to the documentation of this file.
3 
6  : elecsToken_(consumes<std::vector<pat::Electron> >(cfg.getParameter<edm::InputTag>("electrons"))),
7  musToken_(consumes<std::vector<pat::Muon> >(cfg.getParameter<edm::InputTag>("muons"))),
8  jetsToken_(consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
9  metsToken_(consumes<std::vector<pat::MET> >(cfg.getParameter<edm::InputTag>("mets"))),
10 
11  lepton_(nullptr),
12  leptonBar_(nullptr),
13  b_(nullptr),
14  bBar_(nullptr),
15  neutrino_(nullptr),
16  neutrinoBar_(nullptr) {
17  getMatch_ = false;
18  if (cfg.exists("match")) {
19  getMatch_ = true;
20  matchToken_ = consumes<std::vector<std::vector<int> > >(cfg.getParameter<edm::InputTag>("match"));
21  }
22  // if no other correction is given apply L3 (abs) correction
23  jetCorrectionLevel_ = "abs";
24  if (cfg.exists("jetCorrectionLevel")) {
25  jetCorrectionLevel_ = cfg.getParameter<std::string>("jetCorrectionLevel");
26  } else { // if no other correction is given apply L3 (abs) correction
27  jetCorrectionLevel_ = "abs";
28  }
29  produces<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >();
30  produces<int>("Key");
31 }
32 
35  if (lepton_)
36  delete lepton_;
37  if (leptonBar_)
38  delete leptonBar_;
39  if (b_)
40  delete b_;
41  if (bBar_)
42  delete bBar_;
43  if (neutrino_)
44  delete neutrino_;
45  if (neutrinoBar_)
46  delete neutrinoBar_;
47  //if( met_ ) delete met_;
48 }
49 
53  evt.getByToken(elecsToken_, elecs);
54 
56  evt.getByToken(musToken_, mus);
57 
59  evt.getByToken(jetsToken_, jets);
60 
62  evt.getByToken(metsToken_, mets);
63 
64  std::vector<std::vector<int> > matchVec;
65  if (getMatch_) {
67  evt.getByToken(matchToken_, matchHandle);
68  ;
69  matchVec = *matchHandle;
70  } else {
71  std::vector<int> dummyMatch;
72  for (unsigned int i = 0; i < 4; ++i)
73  dummyMatch.push_back(-1);
74  matchVec.push_back(dummyMatch);
75  }
76 
77  // declare unique_ptr for products
78  std::unique_ptr<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > > pOut(
79  new std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > >);
80  std::unique_ptr<int> pKey(new int);
81 
82  // build and feed out key
83  buildKey();
84  *pKey = key();
85  evt.put(std::move(pKey), "Key");
86 
87  // go through given vector of jet combinations
88  unsigned int idMatch = 0;
89  typedef std::vector<std::vector<int> >::iterator MatchVecIterator;
90  for (MatchVecIterator match = matchVec.begin(); match != matchVec.end(); ++match) {
91  // reset pointers
93  // build hypothesis
94  buildHypo(evt, elecs, mus, jets, mets, *match, idMatch++);
95  pOut->push_back(std::make_pair(hypo(), *match));
96  }
97  // feed out hyps and matches
98  evt.put(std::move(pOut));
99 }
100 
103  lepton_ = nullptr;
104  leptonBar_ = nullptr;
105  b_ = nullptr;
106  bBar_ = nullptr;
107  neutrino_ = nullptr;
108  neutrinoBar_ = nullptr;
109  //met_ = 0;
110 }
111 
114  // check for sanity of the hypothesis
115  if (!lepton_ || !leptonBar_ || !b_ || !bBar_) {
116  return reco::CompositeCandidate();
117  }
118 
119  if (key() == TtFullLeptonicEvent::kGenMatch && (!recNu || !recNuBar)) {
120  edm::LogInfo("TtFullHypothesis") << "no neutrinos for gen match" << std::endl;
121  return reco::CompositeCandidate();
122  }
124  edm::LogInfo("TtFullHypothesis") << "no neutrinos for kin solution" << std::endl;
125  return reco::CompositeCandidate();
126  }
127 
128  // setup transient references
130 
131  AddFourMomenta addFourMomenta;
132 
133  // build up the top branch
137  else if (key() == TtFullLeptonicEvent::kGenMatch)
139  addFourMomenta.set(WPlus);
142  addFourMomenta.set(Top);
143 
144  // build up the anti top branch
148  else if (key() == TtFullLeptonicEvent::kGenMatch)
150  addFourMomenta.set(WMinus);
151  TopBar.addDaughter(WMinus, TtFullLepDaughter::WMinus);
153  addFourMomenta.set(TopBar);
154 
155  // build ttbar hypothesis
158  addFourMomenta.set(hyp);
159 
160  // the four momentum of the met is not added to the hypothesis
161  // because it is allready included through the neutrinos
162  //hyp.addDaughter(*met_, TtFullLepDaughter::Met);
163  return hyp;
164 }
165 
167 void TtFullLepHypothesis::setCandidate(const edm::Handle<std::vector<pat::Jet> >& handle,
168  const int& idx,
170  const std::string& correctionLevel) {
172  clone = new reco::ShallowClonePtrCandidate(
173  ptr, ptr->charge(), ptr->correctedJet(jetCorrectionLevel_, "bottom").p4(), ptr->vertex());
174 }
reco::ShallowClonePtrCandidate * neutrino_
reco::ShallowClonePtrCandidate * leptonBar_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
static const std::string WMinus
tuple cfg
Definition: looper.py:296
reco::LeafCandidate * recNuBar
static const std::string Nu
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
virtual void buildKey()=0
build the event hypothesis key
bool exists(std::string const &parameterName) const
checks if a parameter exists
reco::LeafCandidate * recNu
candidates needed for the genmatch hypothesis
static const std::string WPlus
static const std::string B
static const std::string Top
static const std::string Top
static const std::string TopBar
vector< PseudoJet > jets
def move
Definition: eostools.py:511
tuple handle
Definition: patZpeak.py:23
static const std::string WMinus
static const std::string LepBar
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
void resetCandidates()
reset candidate pointers before hypo build process
reco::CompositeCandidate hypo()
return event hypothesis
static const std::string WPlus
Log< level::Info, false > LogInfo
static const std::string TopBar
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_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
reco::ShallowClonePtrCandidate * lepton_
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
reco::ShallowClonePtrCandidate * bBar_
void produce(edm::Event &, const edm::EventSetup &) override
produce the event hypothesis as CompositeCandidate and Key
std::string jetCorrectionLevel_
constexpr char Jet[]
Definition: modules.cc:9
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
static const std::string BBar
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
void set(reco::Candidate &c) const
set up a candidate
constexpr char Electron[]
Definition: modules.cc:12
reco::ShallowClonePtrCandidate * b_
TtFullLepHypothesis(const edm::ParameterSet &)
default constructor
static const std::string NuBar
~TtFullLepHypothesis() override
default destructor
static const std::string Lep
reco::ShallowClonePtrCandidate * neutrinoBar_
edm::EDGetTokenT< std::vector< pat::Muon > > musToken_