CMS 3D CMS Logo

PseudoTopProducer.cc
Go to the documentation of this file.
3 
5 #include "fastjet/ClusterSequence.hh"
6 
11 
14 
15 #include "fastjet/JetDefinition.hh"
16 #include <set>
17 
19 public:
21  void produce(edm::Event& event, const edm::EventSetup& eventSetup) override;
22 
23 private:
24  bool isFromHadron(const reco::Candidate* p) const;
25  bool isBHadron(const reco::Candidate* p) const;
26  bool isBHadron(const unsigned int pdgId) const;
27 
29 
30 private:
34  const double wMass_, tMass_;
40 
41  typedef fastjet::JetDefinition JetDef;
42  std::shared_ptr<JetDef> fjLepDef_, fjJetDef_;
44 };
45 
46 using namespace std;
47 using namespace edm;
48 using namespace reco;
49 
51  : finalStateToken_(consumes<edm::View<reco::Candidate> >(pset.getParameter<edm::InputTag>("finalStates"))),
52  genParticleToken_(consumes<edm::View<reco::Candidate> >(pset.getParameter<edm::InputTag>("genParticles"))),
53  minLeptonPt_(pset.getParameter<double>("minLeptonPt")),
54  maxLeptonEta_(pset.getParameter<double>("maxLeptonEta")),
55  minJetPt_(pset.getParameter<double>("minJetPt")),
56  maxJetEta_(pset.getParameter<double>("maxJetEta")),
57  wMass_(pset.getParameter<double>("wMass")),
58  tMass_(pset.getParameter<double>("tMass")),
59  minLeptonPtDilepton_(pset.getParameter<double>("minLeptonPtDilepton")),
60  maxLeptonEtaDilepton_(pset.getParameter<double>("maxLeptonEtaDilepton")),
61  minDileptonMassDilepton_(pset.getParameter<double>("minDileptonMassDilepton")),
62  minLeptonPtSemilepton_(pset.getParameter<double>("minLeptonPtSemilepton")),
63  maxLeptonEtaSemilepton_(pset.getParameter<double>("maxLeptonEtaSemilepton")),
64  minVetoLeptonPtSemilepton_(pset.getParameter<double>("minVetoLeptonPtSemilepton")),
65  maxVetoLeptonEtaSemilepton_(pset.getParameter<double>("maxVetoLeptonEtaSemilepton")),
66  minMETSemiLepton_(pset.getParameter<double>("minMETSemiLepton")),
67  minMtWSemiLepton_(pset.getParameter<double>("minMtWSemiLepton")) {
68  const double leptonConeSize = pset.getParameter<double>("leptonConeSize");
69  const double jetConeSize = pset.getParameter<double>("jetConeSize");
70  fjLepDef_ = std::make_shared<JetDef>(fastjet::antikt_algorithm, leptonConeSize);
71  fjJetDef_ = std::make_shared<JetDef>(fastjet::antikt_algorithm, jetConeSize);
72 
74 
75  produces<reco::GenParticleCollection>("neutrinos");
76  produces<reco::GenJetCollection>("leptons");
77  produces<reco::GenJetCollection>("jets");
78 
79  produces<reco::GenParticleCollection>();
80 }
81 
83  edm::Handle<edm::View<reco::Candidate> > finalStateHandle;
84  event.getByToken(finalStateToken_, finalStateHandle);
85 
86  edm::Handle<edm::View<reco::Candidate> > genParticleHandle;
87  event.getByToken(genParticleToken_, genParticleHandle);
88 
89  std::unique_ptr<reco::GenParticleCollection> neutrinos(new reco::GenParticleCollection);
90  std::unique_ptr<reco::GenJetCollection> leptons(new reco::GenJetCollection);
91  std::unique_ptr<reco::GenJetCollection> jets(new reco::GenJetCollection);
92  auto neutrinosRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("neutrinos");
93  auto leptonsRefHandle = event.getRefBeforePut<reco::GenJetCollection>("leptons");
94  auto jetsRefHandle = event.getRefBeforePut<reco::GenJetCollection>("jets");
95 
96  std::unique_ptr<reco::GenParticleCollection> pseudoTop(new reco::GenParticleCollection);
97  auto pseudoTopRefHandle = event.getRefBeforePut<reco::GenParticleCollection>();
98 
99  // Collect unstable B-hadrons
100  std::set<int> bHadronIdxs;
101  for (size_t i = 0, n = genParticleHandle->size(); i < n; ++i) {
102  const reco::Candidate& p = genParticleHandle->at(i);
103  const int status = p.status();
104  if (status == 1)
105  continue;
106 
107  // Collect B-hadrons, to be used in b tagging
108  if (isBHadron(&p))
109  bHadronIdxs.insert(-i - 1);
110  }
111 
112  // Collect stable leptons and neutrinos
113  size_t nStables = 0;
114  std::vector<size_t> leptonIdxs;
115  for (size_t i = 0, n = finalStateHandle->size(); i < n; ++i) {
116  const reco::Candidate& p = finalStateHandle->at(i);
117  const int absPdgId = abs(p.pdgId());
118  if (p.status() != 1)
119  continue;
120 
121  ++nStables;
122  if (p.numberOfMothers() == 0)
123  continue; // Skip orphans (if exists)
124  if (p.mother()->status() == 4)
125  continue; // Treat particle as hadronic if directly from the incident beam (protect orphans in MINIAOD)
126  if (isFromHadron(&p))
127  continue;
128  switch (absPdgId) {
129  case 11:
130  case 13: // Leptons
131  case 22: // Photons
132  leptonIdxs.push_back(i);
133  break;
134  case 12:
135  case 14:
136  case 16:
137  neutrinos->push_back(reco::GenParticle(p.charge(), p.p4(), p.vertex(), p.pdgId(), p.status(), true));
138  break;
139  }
140  }
141 
142  // Sort neutrinos by pT.
144 
145  // Make dressed leptons with anti-kt(0.1) algorithm
147  std::vector<fastjet::PseudoJet> fjLepInputs;
148  fjLepInputs.reserve(leptonIdxs.size());
149  for (auto index : leptonIdxs) {
150  const reco::Candidate& p = finalStateHandle->at(index);
151  if (std::isnan(p.pt()) or p.pt() <= 0)
152  continue;
153 
154  fjLepInputs.push_back(fastjet::PseudoJet(p.px(), p.py(), p.pz(), p.energy()));
155  fjLepInputs.back().set_user_index(index);
156  }
157 
159  fastjet::ClusterSequence fjLepClusterSeq(fjLepInputs, *fjLepDef_);
160  std::vector<fastjet::PseudoJet> fjLepJets = fastjet::sorted_by_pt(fjLepClusterSeq.inclusive_jets(minLeptonPt_));
161 
163  leptons->reserve(fjLepJets.size());
164  std::set<size_t> lepDauIdxs; // keep lepton constituents to remove from GenJet construction
165  for (auto& fjJet : fjLepJets) {
166  if (abs(fjJet.eta()) > maxLeptonEta_)
167  continue;
168 
169  // Get jet constituents from fastJet
170  const std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
171  // Convert to CandidatePtr
172  std::vector<reco::CandidatePtr> constituents;
173  reco::CandidatePtr lepCand;
174  for (auto& fjConstituent : fjConstituents) {
175  const size_t index = fjConstituent.user_index();
176  reco::CandidatePtr cand = finalStateHandle->ptrAt(index);
177  const int absPdgId = abs(cand->pdgId());
178  if (absPdgId == 11 or absPdgId == 13) {
179  if (lepCand.isNonnull() and lepCand->pt() > cand->pt())
180  continue; // Choose one with highest pt
181  lepCand = cand;
182  }
183  constituents.push_back(cand);
184  }
185  if (lepCand.isNull())
186  continue;
187 
188  const LorentzVector jetP4(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E());
189  reco::GenJet lepJet;
190  reco::writeSpecific(lepJet, jetP4, genVertex_, constituents);
191 
192  lepJet.setPdgId(lepCand->pdgId());
193  lepJet.setCharge(lepCand->charge());
194 
195  const double jetArea = fjJet.has_area() ? fjJet.area() : 0;
196  lepJet.setJetArea(jetArea);
197 
198  leptons->push_back(lepJet);
199 
200  // Keep constituent indices to be used in the next step.
201  for (auto& fjConstituent : fjConstituents) {
202  lepDauIdxs.insert(fjConstituent.user_index());
203  }
204  }
205 
206  // Now proceed to jets.
207  // Jets: anti-kt excluding the e, mu, nu, and photons in selected leptons.
209  std::vector<fastjet::PseudoJet> fjJetInputs;
210  fjJetInputs.reserve(nStables);
211  for (size_t i = 0, n = finalStateHandle->size(); i < n; ++i) {
212  const reco::Candidate& p = finalStateHandle->at(i);
213  if (p.status() != 1)
214  continue;
215  if (std::isnan(p.pt()) or p.pt() <= 0)
216  continue;
217 
218  const int absId = std::abs(p.pdgId());
219  if (absId == 12 or absId == 14 or absId == 16)
220  continue;
221  if (lepDauIdxs.find(i) != lepDauIdxs.end())
222  continue;
223 
224  fjJetInputs.push_back(fastjet::PseudoJet(p.px(), p.py(), p.pz(), p.energy()));
225  fjJetInputs.back().set_user_index(i);
226  }
228  for (auto index : bHadronIdxs) {
229  const reco::Candidate& p = genParticleHandle->at(abs(index + 1));
230  if (std::isnan(p.pt()) or p.pt() <= 0)
231  continue;
232 
233  const double scale = 1e-20 / p.p();
234  fjJetInputs.push_back(fastjet::PseudoJet(p.px() * scale, p.py() * scale, p.pz() * scale, p.energy() * scale));
235  fjJetInputs.back().set_user_index(index);
236  }
237 
239  fastjet::ClusterSequence fjJetClusterSeq(fjJetInputs, *fjJetDef_);
240  std::vector<fastjet::PseudoJet> fjJets = fastjet::sorted_by_pt(fjJetClusterSeq.inclusive_jets(minJetPt_));
241 
243  jets->reserve(fjJets.size());
244  std::vector<size_t> bjetIdxs, ljetIdxs;
245  for (auto& fjJet : fjJets) {
246  if (abs(fjJet.eta()) > maxJetEta_)
247  continue;
248 
249  // Get jet constituents from fastJet
250  const std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
251  // Convert to CandidatePtr
252  std::vector<reco::CandidatePtr> constituents;
253  bool hasBHadron = false;
254  for (size_t j = 0, m = fjConstituents.size(); j < m; ++j) {
255  const int index = fjConstituents[j].user_index();
256  if (bHadronIdxs.find(index) != bHadronIdxs.end()) {
257  hasBHadron = true;
258  continue;
259  }
260  reco::CandidatePtr cand = finalStateHandle->ptrAt(index);
261  constituents.push_back(cand);
262  }
263 
264  const LorentzVector jetP4(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E());
265  reco::GenJet genJet;
266  reco::writeSpecific(genJet, jetP4, genVertex_, constituents);
267 
268  const double jetArea = fjJet.has_area() ? fjJet.area() : 0;
269  genJet.setJetArea(jetArea);
270  if (hasBHadron) {
271  genJet.setPdgId(5);
272  bjetIdxs.push_back(jets->size());
273  } else {
274  ljetIdxs.push_back(jets->size());
275  }
276 
277  jets->push_back(genJet);
278  }
279 
280  // Every building blocks are ready. Continue to pseudo-W and pseudo-top combination
281  // NOTE : A C++ trick, use do-while instead of long-nested if-statements.
282  do {
283  if (bjetIdxs.size() < 2)
284  break;
285 
286  // Note : we will do dilepton or semilepton channel only
287  if (leptons->size() == 2 and neutrinos->size() >= 2) {
288  // Start from dilepton channel
289  const int q1 = leptons->at(0).charge();
290  const int q2 = leptons->at(1).charge();
291  if (q1 * q2 > 0)
292  break;
293 
294  const auto& lepton1 = q1 > 0 ? leptons->at(0) : leptons->at(1);
295  const auto& lepton2 = q1 > 0 ? leptons->at(1) : leptons->at(0);
296 
297  if (lepton1.pt() < minLeptonPtDilepton_ or std::abs(lepton1.eta()) > maxLeptonEtaDilepton_)
298  break;
299  if (lepton2.pt() < minLeptonPtDilepton_ or std::abs(lepton2.eta()) > maxLeptonEtaDilepton_)
300  break;
301  if ((lepton1.p4() + lepton2.p4()).mass() < minDileptonMassDilepton_)
302  break;
303 
304  double dm = 1e9;
305  int selNu1 = -1, selNu2 = -1;
306  for (int i = 0; i < 2; ++i) { // Consider leading 2 neutrinos. Neutrino vector is already sorted by pT
307  const double dm1 = std::abs((lepton1.p4() + neutrinos->at(i).p4()).mass() - wMass_);
308  for (int j = 0; j < 2; ++j) {
309  if (i == j)
310  continue;
311  const double dm2 = std::abs((lepton2.p4() + neutrinos->at(j).p4()).mass() - wMass_);
312  const double newDm = dm1 + dm2;
313 
314  if (newDm < dm) {
315  dm = newDm;
316  selNu1 = i;
317  selNu2 = j;
318  }
319  }
320  }
321  if (dm >= 1e9)
322  break;
323 
324  const auto& nu1 = neutrinos->at(selNu1);
325  const auto& nu2 = neutrinos->at(selNu2);
326  const auto w1LVec = lepton1.p4() + nu1.p4();
327  const auto w2LVec = lepton2.p4() + nu2.p4();
328 
329  // Contiue to top quarks
330  dm = 1e9; // Reset once again for top combination.
331  int selB1 = -1, selB2 = -1;
332  for (auto i : bjetIdxs) {
333  const double dm1 = std::abs((w1LVec + jets->at(i).p4()).mass() - tMass_);
334  for (auto j : bjetIdxs) {
335  if (i == j)
336  continue;
337  const double dm2 = std::abs((w2LVec + jets->at(j).p4()).mass() - tMass_);
338  const double newDm = dm1 + dm2;
339 
340  if (newDm < dm) {
341  dm = newDm;
342  selB1 = i;
343  selB2 = j;
344  }
345  }
346  }
347  if (dm >= 1e9)
348  break;
349 
350  const auto& bJet1 = jets->at(selB1);
351  const auto& bJet2 = jets->at(selB2);
352  const auto t1LVec = w1LVec + bJet1.p4();
353  const auto t2LVec = w2LVec + bJet2.p4();
354 
355  // Recalculate lepton charges (needed due to initialization of leptons wrt sign of a charge)
356  const int lepQ1 = lepton1.charge();
357  const int lepQ2 = lepton2.charge();
358 
359  // Put all of them into candidate collection
360  reco::GenParticle t1(lepQ1 * 2 / 3., t1LVec, genVertex_, lepQ1 * 6, 3, false);
361  reco::GenParticle w1(lepQ1, w1LVec, genVertex_, lepQ1 * 24, 3, true);
362  reco::GenParticle b1(-lepQ1 / 3., bJet1.p4(), genVertex_, lepQ1 * 5, 1, true);
363  reco::GenParticle l1(lepQ1, lepton1.p4(), genVertex_, lepton1.pdgId(), 1, true);
364  reco::GenParticle n1(0, nu1.p4(), genVertex_, nu1.pdgId(), 1, true);
365 
366  reco::GenParticle t2(lepQ2 * 2 / 3., t2LVec, genVertex_, lepQ2 * 6, 3, false);
367  reco::GenParticle w2(lepQ2, w2LVec, genVertex_, lepQ2 * 24, 3, true);
368  reco::GenParticle b2(-lepQ2 / 3., bJet2.p4(), genVertex_, lepQ2 * 5, 1, true);
369  reco::GenParticle l2(lepQ2, lepton2.p4(), genVertex_, lepton2.pdgId(), 1, true);
370  reco::GenParticle n2(0, nu2.p4(), genVertex_, nu2.pdgId(), 1, true);
371 
372  pseudoTop->push_back(t1);
373  pseudoTop->push_back(t2);
374 
375  pseudoTop->push_back(w1);
376  pseudoTop->push_back(b1);
377  pseudoTop->push_back(l1);
378  pseudoTop->push_back(n1);
379 
380  pseudoTop->push_back(w2);
381  pseudoTop->push_back(b2);
382  pseudoTop->push_back(l2);
383  pseudoTop->push_back(n2);
384  } else if (!leptons->empty() and !neutrinos->empty() and ljetIdxs.size() >= 2) {
385  // Then continue to the semi-leptonic channel
386  const auto& lepton = leptons->at(0);
387  if (lepton.pt() < minLeptonPtSemilepton_ or std::abs(lepton.eta()) > maxLeptonEtaSemilepton_)
388  break;
389 
390  // Skip event if there are veto leptons
391  bool hasVetoLepton = false;
392  for (auto vetoLeptonItr = next(leptons->begin()); vetoLeptonItr != leptons->end(); ++vetoLeptonItr) {
393  if (vetoLeptonItr->pt() > minVetoLeptonPtSemilepton_ and
394  std::abs(vetoLeptonItr->eta()) < maxVetoLeptonEtaSemilepton_) {
395  hasVetoLepton = true;
396  }
397  }
398  if (hasVetoLepton)
399  break;
400 
401  // Calculate MET
402  double metX = 0, metY = 0;
403  for (const auto& neutrino : *neutrinos) {
404  metX += neutrino.px();
405  metY += neutrino.py();
406  }
407  const double metPt = std::hypot(metX, metY);
408  if (metPt < minMETSemiLepton_)
409  break;
410 
411  const double mtW = std::sqrt(2 * lepton.pt() * metPt * cos(lepton.phi() - atan2(metX, metY)));
412  if (mtW < minMtWSemiLepton_)
413  break;
414 
415  // Solve pz to give wMass_^2 = (lepE+energy)^2 - (lepPx+metX)^2 - (lepPy+metY)^2 - (lepPz+pz)^2
416  // -> (wMass_^2)/2 = lepE*sqrt(metPt^2+pz^2) - lepPx*metX - lepPy*metY - lepPz*pz
417  // -> C := (wMass_^2)/2 + lepPx*metX + lepPy*metY
418  // -> lepE^2*(metPt^2+pz^2) = C^2 + 2*C*lepPz*pz + lepPz^2*pz^2
419  // -> (lepE^2-lepPz^2)*pz^2 - 2*C*lepPz*pz - C^2 + lepE^2*metPt^2 = 0
420  // -> lepPt^2*pz^2 - 2*C*lepPz*pz - C^2 + lepE^2*metPt^2 = 0
421  const double lpz = lepton.pz(), le = lepton.energy(), lpt = lepton.pt();
422  const double cf = (wMass_ * wMass_) / 2 + lepton.px() * metX + lepton.py() * metY;
423  const double det = cf * cf * lpz * lpz - lpt * lpt * (le * le * metPt * metPt - cf * cf);
424  const double pz = (cf * lpz + (cf < 0 ? -sqrt(det) : sqrt(det))) / lpt / lpt;
425  const reco::Candidate::LorentzVector nu1P4(metX, metY, pz, std::hypot(metPt, pz));
426  const auto w1LVec = lepton.p4() + nu1P4;
427 
428  // Continue to build leptonic pseudo top
429  double minDR = 1e9;
430  int selB1 = -1;
431  for (auto b1Itr = bjetIdxs.begin(); b1Itr != bjetIdxs.end(); ++b1Itr) {
432  const double dR = deltaR(jets->at(*b1Itr).p4(), w1LVec);
433  if (dR < minDR) {
434  selB1 = *b1Itr;
435  minDR = dR;
436  }
437  }
438  if (selB1 == -1)
439  break;
440  const auto& bJet1 = jets->at(selB1);
441  const auto t1LVec = w1LVec + bJet1.p4();
442 
443  // Build hadronic pseudo W, take leading two jets (ljetIdxs are (implicitly) sorted by pT)
444  const auto& wJet1 = jets->at(ljetIdxs[0]);
445  const auto& wJet2 = jets->at(ljetIdxs[1]);
446  const auto w2LVec = wJet1.p4() + wJet2.p4();
447 
448  // Contiue to top quarks
449  double dm = 1e9;
450  int selB2 = -1;
451  for (auto i : bjetIdxs) {
452  if (int(i) == selB1)
453  continue;
454  const double newDm = std::abs((w2LVec + jets->at(i).p4()).mass() - tMass_);
455  if (newDm < dm) {
456  dm = newDm;
457  selB2 = i;
458  }
459  }
460  if (dm >= 1e9)
461  break;
462 
463  const auto& bJet2 = jets->at(selB2);
464  const auto t2LVec = w2LVec + bJet2.p4();
465 
466  const int q = lepton.charge();
467  // Put all of them into candidate collection
468  reco::GenParticle t1(q * 2 / 3., t1LVec, genVertex_, q * 6, 3, false);
469  reco::GenParticle w1(q, w1LVec, genVertex_, q * 24, 3, true);
470  reco::GenParticle b1(-q / 3., bJet1.p4(), genVertex_, q * 5, 1, true);
471  reco::GenParticle l1(q, lepton.p4(), genVertex_, lepton.pdgId(), 1, true);
472  reco::GenParticle n1(0, nu1P4, genVertex_, q * (std::abs(lepton.pdgId()) + 1), 1, true);
473 
474  reco::GenParticle t2(-q * 2 / 3., t2LVec, genVertex_, -q * 6, 3, false);
475  reco::GenParticle w2(-q, w2LVec, genVertex_, -q * 24, 3, true);
476  reco::GenParticle b2(0, bJet2.p4(), genVertex_, -q * 5, 1, true);
477  reco::GenParticle u2(0, wJet1.p4(), genVertex_, -2 * q, 1, true);
478  reco::GenParticle d2(0, wJet2.p4(), genVertex_, q, 1, true);
479 
480  pseudoTop->push_back(t1);
481  pseudoTop->push_back(t2);
482 
483  pseudoTop->push_back(w1);
484  pseudoTop->push_back(b1);
485  pseudoTop->push_back(l1);
486  pseudoTop->push_back(n1);
487 
488  pseudoTop->push_back(w2);
489  pseudoTop->push_back(b2);
490  pseudoTop->push_back(u2);
491  pseudoTop->push_back(d2);
492  }
493  } while (false);
494 
495  if (pseudoTop->size() == 10) // If pseudtop decay tree is completed
496  {
497  // t->W+b
498  pseudoTop->at(0).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 2)); // t->W
499  pseudoTop->at(0).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 3)); // t->b
500  pseudoTop->at(2).addMother(reco::GenParticleRef(pseudoTopRefHandle, 0)); // t->W
501  pseudoTop->at(3).addMother(reco::GenParticleRef(pseudoTopRefHandle, 0)); // t->b
502 
503  // W->lv or W->jj
504  pseudoTop->at(2).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 4));
505  pseudoTop->at(2).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 5));
506  pseudoTop->at(4).addMother(reco::GenParticleRef(pseudoTopRefHandle, 2));
507  pseudoTop->at(5).addMother(reco::GenParticleRef(pseudoTopRefHandle, 2));
508 
509  // tbar->W-b
510  pseudoTop->at(1).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 6));
511  pseudoTop->at(1).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 7));
512  pseudoTop->at(6).addMother(reco::GenParticleRef(pseudoTopRefHandle, 1));
513  pseudoTop->at(7).addMother(reco::GenParticleRef(pseudoTopRefHandle, 1));
514 
515  // W->jj
516  pseudoTop->at(6).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 8));
517  pseudoTop->at(6).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 9));
518  pseudoTop->at(8).addMother(reco::GenParticleRef(pseudoTopRefHandle, 6));
519  pseudoTop->at(9).addMother(reco::GenParticleRef(pseudoTopRefHandle, 6));
520  }
521 
522  event.put(std::move(neutrinos), "neutrinos");
523  event.put(std::move(leptons), "leptons");
524  event.put(std::move(jets), "jets");
525 
526  event.put(std::move(pseudoTop));
527 }
528 
530  for (size_t i = 0, n = p->numberOfMothers(); i < n; ++i) {
531  const reco::Candidate* mother = p->mother(i);
532  if (!mother)
533  continue;
534  if (mother->status() == 4 or mother->numberOfMothers() == 0)
535  continue; // Skip incident beam
536  const int pdgId = abs(mother->pdgId());
537 
538  if (pdgId > 100)
539  return true;
540  else if (isFromHadron(mother))
541  return true;
542  }
543  return false;
544 }
545 
547  const unsigned int absPdgId = abs(p->pdgId());
548  if (!isBHadron(absPdgId))
549  return false;
550 
551  // Do not consider this particle if it has B hadron daughter
552  // For example, B* -> B0 + photon; then we drop B* and take B0 only
553  for (int i = 0, n = p->numberOfDaughters(); i < n; ++i) {
554  const reco::Candidate* dau = p->daughter(i);
555  if (isBHadron(abs(dau->pdgId())))
556  return false;
557  }
558 
559  return true;
560 }
561 
562 bool PseudoTopProducer::isBHadron(const unsigned int absPdgId) const {
563  if (absPdgId <= 100)
564  return false; // Fundamental particles and MC internals
565  if (absPdgId >= 1000000000)
566  return false; // Nuclei, +-10LZZZAAAI
567 
568  // General form of PDG ID is 7 digit form
569  // +- n nr nL nq1 nq2 nq3 nJ
570  //const int nJ = absPdgId % 10; // Spin
571  const int nq3 = (absPdgId / 10) % 10;
572  const int nq2 = (absPdgId / 100) % 10;
573  const int nq1 = (absPdgId / 1000) % 10;
574 
575  if (nq3 == 0)
576  return false; // Diquarks
577  if (nq1 == 0 and nq2 == 5)
578  return true; // B mesons
579  if (nq1 == 5)
580  return true; // B baryons
581 
582  return false;
583 }
584 
const edm::EDGetTokenT< edm::View< reco::Candidate > > genParticleToken_
PseudoTopProducer(const edm::ParameterSet &pset)
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
def isnan(num)
reco::Particle::Point genVertex_
MPlex< T, D1, D2, N > hypot(const MPlex< T, D1, D2, N > &a, const MPlex< T, D1, D2, N > &b)
Definition: Matriplex.h:436
const double maxLeptonEtaSemilepton_
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
const double maxLeptonEta_
const double minVetoLeptonPtSemilepton_
fastjet::JetDefinition JetDef
const double maxLeptonEtaDilepton_
virtual double pt() const =0
transverse momentum
std::shared_ptr< JetDef > fjJetDef_
virtual int status() const =0
status word
std::vector< GenJet > GenJetCollection
collection of GenJet objects
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, CaloGeometry const &geometry, HcalTopology const &topology)
Definition: JetSpecific.cc:32
std::shared_ptr< JetDef > fjLepDef_
const double minMtWSemiLepton_
const double minMETSemiLepton_
reco::Particle::LorentzVector LorentzVector
const double maxJetEta_
const double minLeptonPtSemilepton_
bool isNull() const
Checks for null.
Definition: Ptr.h:142
bJet1
do the razor with one or two b jets (medium CSV)
T sqrt(T t)
Definition: SSEVec.h:19
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:146
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
math::XYZPoint Point
point in the space
Definition: Particle.h:25
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:23
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual int charge() const =0
electric charge
bias2_t b2[25]
Definition: b2.h:9
virtual int pdgId() const =0
PDG identifier.
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
bool isFromHadron(const reco::Candidate *p) const
fixed size matrix
HLT enums.
const double maxVetoLeptonEtaSemilepton_
const double minLeptonPt_
const double minDileptonMassDilepton_
const edm::EDGetTokenT< edm::View< reco::Candidate > > finalStateToken_
const double minLeptonPtDilepton_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
const double minJetPt_
def move(src, dest)
Definition: eostools.py:511
static constexpr float b1
Definition: event.py:1
bool isBHadron(const reco::Candidate *p) const