test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
PseudoTopProducer Class Reference

#include <PseudoTopProducer.h>

Inheritance diagram for PseudoTopProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void produce (edm::Event &event, const edm::EventSetup &eventSetup) override
 
 PseudoTopProducer (const edm::ParameterSet &pset)
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (std::string const &iProcessName, std::string const &iModuleLabel, bool iPrint, std::vector< char const * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Types

typedef fastjet::JetDefinition JetDef
 
typedef
reco::Particle::LorentzVector 
LorentzVector
 

Private Member Functions

reco::GenParticleRef buildGenParticle (const reco::Candidate *p, reco::GenParticleRefProd &refHandle, std::auto_ptr< reco::GenParticleCollection > &outColl) const
 
const reco::CandidategetLast (const reco::Candidate *p)
 
bool isBHadron (const reco::Candidate *p) const
 
bool isBHadron (const unsigned int pdgId) const
 
bool isFromHadron (const reco::Candidate *p) const
 

Private Attributes

edm::EDGetTokenT< edm::View
< reco::Candidate > > 
finalStateToken_
 
std::shared_ptr< JetDeffjJetDef_
 
std::shared_ptr< JetDeffjLepDef_
 
edm::EDGetTokenT< edm::View
< reco::Candidate > > 
genParticleToken_
 
reco::Particle::Point genVertex_
 
const double jetMaxEta_
 
const double jetMinPt_
 
const double leptonMaxEta_
 
const double leptonMinPt_
 
const double tMass_
 
const double wMass_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 13 of file PseudoTopProducer.h.

Member Typedef Documentation

typedef fastjet::JetDefinition PseudoTopProducer::JetDef
private

Definition at line 36 of file PseudoTopProducer.h.

Definition at line 28 of file PseudoTopProducer.h.

Constructor & Destructor Documentation

PseudoTopProducer::PseudoTopProducer ( const edm::ParameterSet pset)

Definition at line 12 of file PseudoTopProducer.cc.

References finalStateToken_, fjJetDef_, fjLepDef_, genParticleToken_, genVertex_, and edm::ParameterSet::getParameter().

12  :
13  leptonMinPt_(pset.getParameter<double>("leptonMinPt")),
14  leptonMaxEta_(pset.getParameter<double>("leptonMaxEta")),
15  jetMinPt_(pset.getParameter<double>("jetMinPt")),
16  jetMaxEta_(pset.getParameter<double>("jetMaxEta")),
17  wMass_(pset.getParameter<double>("wMass")),
18  tMass_(pset.getParameter<double>("tMass"))
19 {
20  finalStateToken_ = consumes<edm::View<reco::Candidate> >(pset.getParameter<edm::InputTag>("finalStates"));
21  genParticleToken_ = consumes<edm::View<reco::Candidate> >(pset.getParameter<edm::InputTag>("genParticles"));
22 
23  const double leptonConeSize = pset.getParameter<double>("leptonConeSize");
24  const double jetConeSize = pset.getParameter<double>("jetConeSize");
25  fjLepDef_ = std::shared_ptr<JetDef>(new JetDef(fastjet::antikt_algorithm, leptonConeSize));
26  fjJetDef_ = std::shared_ptr<JetDef>(new JetDef(fastjet::antikt_algorithm, jetConeSize));
27 
29 
30  produces<reco::GenParticleCollection>("neutrinos");
31  produces<reco::GenJetCollection>("leptons");
32  produces<reco::GenJetCollection>("jets");
33 
34  produces<reco::GenParticleCollection>();
35 
36 }
const double leptonMaxEta_
T getParameter(std::string const &) const
reco::Particle::Point genVertex_
const double leptonMinPt_
const double jetMaxEta_
edm::EDGetTokenT< edm::View< reco::Candidate > > finalStateToken_
fastjet::JetDefinition JetDef
std::shared_ptr< JetDef > fjJetDef_
std::shared_ptr< JetDef > fjLepDef_
math::XYZPoint Point
point in the space
Definition: Particle.h:25
edm::EDGetTokenT< edm::View< reco::Candidate > > genParticleToken_
const double jetMinPt_

Member Function Documentation

reco::GenParticleRef PseudoTopProducer::buildGenParticle ( const reco::Candidate p,
reco::GenParticleRefProd refHandle,
std::auto_ptr< reco::GenParticleCollection > &  outColl 
) const
private

Definition at line 566 of file PseudoTopProducer.cc.

References reco::CompositeRefCandidateT< D >::clearDaughters(), reco::CompositeRefCandidateT< D >::clearMothers(), edm::RefProd< T >::id(), reco::CompositeRefCandidateT< D >::resetDaughters(), and reco::CompositeRefCandidateT< D >::resetMothers().

568 {
569  reco::GenParticle pOut(*dynamic_cast<const reco::GenParticle*>(p));
570  pOut.clearMothers();
571  pOut.clearDaughters();
572  pOut.resetMothers(refHandle.id());
573  pOut.resetDaughters(refHandle.id());
574 
575  outColl->push_back(pOut);
576 
577  return reco::GenParticleRef(refHandle, outColl->size()-1);
578 }
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:137
const reco::Candidate * PseudoTopProducer::getLast ( const reco::Candidate p)
private

Definition at line 507 of file PseudoTopProducer.cc.

References reco::Candidate::daughter(), i, gen::n, reco::Candidate::numberOfDaughters(), AlCaHLTBitMon_ParallelJobs::p, and reco::Candidate::pdgId().

508 {
509  for ( size_t i=0, n=p->numberOfDaughters(); i<n; ++i )
510  {
511  const reco::Candidate* dau = p->daughter(i);
512  if ( p->pdgId() == dau->pdgId() ) return getLast(dau);
513  }
514  return p;
515 }
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual size_type numberOfDaughters() const =0
number of daughters
const reco::Candidate * getLast(const reco::Candidate *p)
virtual int pdgId() const =0
PDG identifier.
bool PseudoTopProducer::isBHadron ( const reco::Candidate p) const
private

Definition at line 531 of file PseudoTopProducer.cc.

References funct::abs(), reco::Candidate::daughter(), i, gen::n, reco::Candidate::numberOfDaughters(), and reco::Candidate::pdgId().

Referenced by produce().

532 {
533  const unsigned int absPdgId = abs(p->pdgId());
534  if ( !isBHadron(absPdgId) ) return false;
535 
536  // Do not consider this particle if it has B hadron daughter
537  // For example, B* -> B0 + photon; then we drop B* and take B0 only
538  for ( int i=0, n=p->numberOfDaughters(); i<n; ++i )
539  {
540  const reco::Candidate* dau = p->daughter(i);
541  if ( isBHadron(abs(dau->pdgId())) ) return false;
542  }
543 
544  return true;
545 }
bool isBHadron(const reco::Candidate *p) const
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual size_type numberOfDaughters() const =0
number of daughters
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int pdgId() const =0
PDG identifier.
bool PseudoTopProducer::isBHadron ( const unsigned int  pdgId) const
private

Definition at line 547 of file PseudoTopProducer.cc.

548 {
549  if ( absPdgId <= 100 ) return false; // Fundamental particles and MC internals
550  if ( absPdgId >= 1000000000 ) return false; // Nuclei, +-10LZZZAAAI
551 
552  // General form of PDG ID is 7 digit form
553  // +- n nr nL nq1 nq2 nq3 nJ
554  //const int nJ = absPdgId % 10; // Spin
555  const int nq3 = (absPdgId / 10) % 10;
556  const int nq2 = (absPdgId / 100) % 10;
557  const int nq1 = (absPdgId / 1000) % 10;
558 
559  if ( nq3 == 0 ) return false; // Diquarks
560  if ( nq1 == 0 and nq2 == 5 ) return true; // B mesons
561  if ( nq1 == 5 ) return true; // B baryons
562 
563  return false;
564 }
bool PseudoTopProducer::isFromHadron ( const reco::Candidate p) const
private

Definition at line 517 of file PseudoTopProducer.cc.

References funct::abs(), i, reco::Candidate::mother(), gen::n, reco::Candidate::numberOfMothers(), SingleMuPt40Fwdv2_cfi_GEN_SIM::pdgId, and reco::Candidate::pdgId().

Referenced by produce().

518 {
519  for ( size_t i=0, n=p->numberOfMothers(); i<n; ++i )
520  {
521  const reco::Candidate* mother = p->mother(i);
522  if ( mother->numberOfMothers() == 0 ) continue; // Skip incident beam
523  const int pdgId = abs(mother->pdgId());
524 
525  if ( pdgId > 100 ) return true;
526  else if ( isFromHadron(mother) ) return true;
527  }
528  return false;
529 }
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
bool isFromHadron(const reco::Candidate *p) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int pdgId() const =0
PDG identifier.
void PseudoTopProducer::produce ( edm::Event event,
const edm::EventSetup eventSetup 
)
overridevirtual

Build jets

Implements edm::EDProducer.

Definition at line 38 of file PseudoTopProducer.cc.

References funct::abs(), reco::CompositeRefCandidateT< D >::addDaughter(), RazorAnalyzer::bJet1, RazorAnalyzer::bJet2, reco::Candidate::charge(), symbols::dm, alignCSCRings::e, reco::Candidate::energy(), finalStateToken_, fjJetDef_, fjLepDef_, genParticleToken_, genVertex_, i, cmsHarvester::index, isBHadron(), isFromHadron(), edm::detail::isnan(), edm::Ptr< T >::isNonnull(), edm::Ptr< T >::isNull(), j, jetMaxEta_, jetMinPt_, fwrapper::jets, leptonMaxEta_, leptonMinPt_, HLT_25ns10e33_v2_cff::leptons, visualization-live-secondInstance_cfg::m, ResonanceBuilder::mass, reco::Candidate::mother(), gen::n, reco::Candidate::numberOfMothers(), or, reco::Candidate::p(), AlCaHLTBitMon_ParallelJobs::p, reco::Candidate::p4(), reco::Candidate::pdgId(), pseudoTop_cfi::pseudoTop, reco::Candidate::pt(), reco::Candidate::px(), reco::Candidate::py(), reco::Candidate::pz(), pileupReCalc_HLTpaths::scale, mps_update::status, reco::Candidate::status(), reco::t2, tMass_, reco::Candidate::vertex(), w2, wMass_, and reco::writeSpecific().

39 {
40  edm::Handle<edm::View<reco::Candidate> > finalStateHandle;
41  event.getByToken(finalStateToken_, finalStateHandle);
42 
43  edm::Handle<edm::View<reco::Candidate> > genParticleHandle;
44  event.getByToken(genParticleToken_, genParticleHandle);
45 
46  std::auto_ptr<reco::GenParticleCollection> neutrinos(new reco::GenParticleCollection);
47  std::auto_ptr<reco::GenJetCollection> leptons(new reco::GenJetCollection);
48  std::auto_ptr<reco::GenJetCollection> jets(new reco::GenJetCollection);
49  auto neutrinosRefHandle = event.getRefBeforePut<reco::GenParticleCollection>("neutrinos");
50  auto leptonsRefHandle = event.getRefBeforePut<reco::GenJetCollection>("leptons");
51  auto jetsRefHandle = event.getRefBeforePut<reco::GenJetCollection>("jets");
52 
53  std::auto_ptr<reco::GenParticleCollection> pseudoTop(new reco::GenParticleCollection);
54  auto pseudoTopRefHandle = event.getRefBeforePut<reco::GenParticleCollection>();
55 
56  // Collect unstable B-hadrons
57  std::set<int> bHadronIdxs;
58  for ( size_t i=0, n=genParticleHandle->size(); i<n; ++i )
59  {
60  const reco::Candidate& p = genParticleHandle->at(i);
61  const int status = p.status();
62  if ( status == 1 ) continue;
63 
64  // Collect B-hadrons, to be used in b tagging
65  if ( isBHadron(&p) ) bHadronIdxs.insert(-i-1);
66  }
67 
68  // Collect stable leptons and neutrinos
69  size_t nStables = 0;
70  std::vector<size_t> leptonIdxs;
71  std::set<size_t> neutrinoIdxs;
72  for ( size_t i=0, n=finalStateHandle->size(); i<n; ++i )
73  {
74  const reco::Candidate& p = finalStateHandle->at(i);
75  const int absPdgId = abs(p.pdgId());
76  if ( p.status() != 1 ) continue;
77 
78  ++nStables;
79  if ( p.numberOfMothers() == 0 ) continue; // Skip orphans (if exists)
80  if ( p.mother()->status() == 4 ) continue; // Treat particle as hadronic if directly from the incident beam (protect orphans in MINIAOD)
81  if ( isFromHadron(&p) ) continue;
82  switch ( absPdgId )
83  {
84  case 11: case 13: // Leptons
85  case 22: // Photons
86  leptonIdxs.push_back(i);
87  break;
88  case 12: case 14: case 16:
89  neutrinoIdxs.insert(i);
90  neutrinos->push_back(reco::GenParticle(p.charge(), p.p4(), p.vertex(), p.pdgId(), p.status(), true));
91  break;
92  }
93  }
94 
95  // Sort neutrinos by pT.
96  std::sort(neutrinos->begin(), neutrinos->end(), GreaterByPt<reco::Candidate>());
97 
98  // Make dressed leptons with anti-kt(0.1) algorithm
100  std::vector<fastjet::PseudoJet> fjLepInputs;
101  fjLepInputs.reserve(leptonIdxs.size());
102  for ( auto index : leptonIdxs )
103  {
104  const reco::Candidate& p = finalStateHandle->at(index);
105  if ( std::isnan(p.pt()) or p.pt() <= 0 ) continue;
106 
107  fjLepInputs.push_back(fastjet::PseudoJet(p.px(), p.py(), p.pz(), p.energy()));
108  fjLepInputs.back().set_user_index(index);
109  }
110 
112  fastjet::ClusterSequence fjLepClusterSeq(fjLepInputs, *fjLepDef_);
113  std::vector<fastjet::PseudoJet> fjLepJets = fastjet::sorted_by_pt(fjLepClusterSeq.inclusive_jets(leptonMinPt_));
114 
116  leptons->reserve(fjLepJets.size());
117  std::set<size_t> lepDauIdxs; // keep lepton constituents to remove from GenJet construction
118  for ( auto& fjJet : fjLepJets )
119  {
120  if ( abs(fjJet.eta()) > leptonMaxEta_ ) continue;
121 
122  // Get jet constituents from fastJet
123  const std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
124  // Convert to CandidatePtr
125  std::vector<reco::CandidatePtr> constituents;
126  reco::CandidatePtr lepCand;
127  for ( auto& fjConstituent : fjConstituents )
128  {
129  const size_t index = fjConstituent.user_index();
130  reco::CandidatePtr cand = finalStateHandle->ptrAt(index);
131  const int absPdgId = abs(cand->pdgId());
132  if ( absPdgId == 11 or absPdgId == 13 )
133  {
134  if ( lepCand.isNonnull() and lepCand->pt() > cand->pt() ) continue; // Choose one with highest pt
135  lepCand = cand;
136  }
137  constituents.push_back(cand);
138  }
139  if ( lepCand.isNull() ) continue;
140  if ( lepCand->pt() < fjJet.pt()/2 ) continue; // Central lepton must be the major component
141 
142  const LorentzVector jetP4(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E());
143  reco::GenJet lepJet;
144  reco::writeSpecific(lepJet, jetP4, genVertex_, constituents, eventSetup);
145 
146  lepJet.setPdgId(lepCand->pdgId());
147  lepJet.setCharge(lepCand->charge());
148 
149  const double jetArea = fjJet.has_area() ? fjJet.area() : 0;
150  lepJet.setJetArea(jetArea);
151 
152  leptons->push_back(lepJet);
153 
154  // Keep constituent indices to be used in the next step.
155  for ( auto& fjConstituent : fjConstituents )
156  {
157  lepDauIdxs.insert(fjConstituent.user_index());
158  }
159  }
160 
161  // Now proceed to jets.
162  // Jets: anti-kt excluding the e, mu, nu, and photons in selected leptons.
164  std::vector<fastjet::PseudoJet> fjJetInputs;
165  fjJetInputs.reserve(nStables);
166  for ( size_t i=0, n=finalStateHandle->size(); i<n; ++i )
167  {
168  const reco::Candidate& p = finalStateHandle->at(i);
169  if ( p.status() != 1 ) continue;
170  if ( std::isnan(p.pt()) or p.pt() <= 0 ) continue;
171 
172  if ( neutrinoIdxs.find(i) != neutrinoIdxs.end() ) continue;
173  if ( lepDauIdxs.find(i) != lepDauIdxs.end() ) continue;
174 
175  fjJetInputs.push_back(fastjet::PseudoJet(p.px(), p.py(), p.pz(), p.energy()));
176  fjJetInputs.back().set_user_index(i);
177  }
179  for ( auto index : bHadronIdxs )
180  {
181  const reco::Candidate& p = genParticleHandle->at(abs(index+1));
182  if ( std::isnan(p.pt()) or p.pt() <= 0 ) continue;
183 
184  const double scale = 1e-20/p.p();
185  fjJetInputs.push_back(fastjet::PseudoJet(p.px()*scale, p.py()*scale, p.pz()*scale, p.energy()*scale));
186  fjJetInputs.back().set_user_index(index);
187  }
188 
190  fastjet::ClusterSequence fjJetClusterSeq(fjJetInputs, *fjJetDef_);
191  std::vector<fastjet::PseudoJet> fjJets = fastjet::sorted_by_pt(fjJetClusterSeq.inclusive_jets(jetMinPt_));
192 
194  jets->reserve(fjJets.size());
195  std::vector<size_t> bjetIdxs, ljetIdxs;
196  for ( auto& fjJet : fjJets )
197  {
198  if ( abs(fjJet.eta()) > jetMaxEta_ ) continue;
199 
200  // Get jet constituents from fastJet
201  const std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(fjJet.constituents());
202  // Convert to CandidatePtr
203  std::vector<reco::CandidatePtr> constituents;
204  bool hasBHadron = false;
205  for ( size_t j=0, m=fjConstituents.size(); j<m; ++j )
206  {
207  const int index = fjConstituents[j].user_index();
208  if ( bHadronIdxs.find(index) != bHadronIdxs.end() ) {
209  hasBHadron = true;
210  continue;
211  }
212  reco::CandidatePtr cand = finalStateHandle->ptrAt(index);
213  constituents.push_back(cand);
214  }
215 
216  const LorentzVector jetP4(fjJet.px(), fjJet.py(), fjJet.pz(), fjJet.E());
217  reco::GenJet genJet;
218  reco::writeSpecific(genJet, jetP4, genVertex_, constituents, eventSetup);
219 
220  const double jetArea = fjJet.has_area() ? fjJet.area() : 0;
221  genJet.setJetArea(jetArea);
222  if ( hasBHadron )
223  {
224  genJet.setPdgId(5);
225  bjetIdxs.push_back(jets->size());
226  }
227  else
228  {
229  ljetIdxs.push_back(jets->size());
230  }
231 
232  jets->push_back(genJet);
233  }
234 
235  // Every building blocks are ready. Continue to pseudo-W and pseudo-top combination
236  // NOTE : A C++ trick, use do-while instead of long-nested if-statements.
237  do
238  {
239  // Note : we will do dilepton or semilepton channel only
240  const size_t nLepton = leptons->size();
241 
242  // Collect leptonic-decaying W's
243  std::vector<std::pair<int, int> > wCandIdxs;
244  for ( auto lep = leptons->begin(); lep != leptons->end(); ++lep )
245  {
246  const size_t iLep = lep-leptons->begin();
247  for ( auto nu = neutrinos->begin(); nu != neutrinos->end(); ++nu )
248  {
249  //if ( abs(lep->pdgId()+nu->pdgId()) != 1 ) continue; // Enforce to conserve flavour, this reduces combinatorial bkg
250  const double m = (lep->p4()+nu->p4()).mass();
251  if ( m > 300 ) continue; // Raw mass cut, reduce combinatorial background
252 
253  const size_t iNu = nu-neutrinos->begin();
254  wCandIdxs.push_back(make_pair(iLep, iNu));
255  }
256  }
257 
258  if ( nLepton == 0 or wCandIdxs.empty() ) break; // Skip full hadronic channel
259  else if ( nLepton == 1 and wCandIdxs.size() >= 1 ) // Semi-leptonic channel
260  {
261  double dm = 1e9; // Note: this will be re-used later.
262  int bestLepIdx = -1, bestNuIdx = -1;
263  for ( auto wCandIdx : wCandIdxs )
264  {
265  const int lepIdx = wCandIdx.first;
266  const int nuIdx = wCandIdx.second;
267  const LorentzVector& lepLVec = leptons->at(lepIdx).p4();
268  const LorentzVector& nuLVec = neutrinos->at(nuIdx).p4();
269  const double mW = (lepLVec + nuLVec).mass();
270  if ( mW > 300 ) continue;
271 
272  const double dmNew = abs(mW-wMass_);
273  if ( dmNew < dm )
274  {
275  dm = dmNew;
276  bestLepIdx = lepIdx;
277  bestNuIdx = nuIdx;
278  }
279  }
280  if ( bestLepIdx < 0 ) break; // Actually this should never happen
281  const auto& lepton = leptons->at(bestLepIdx);
282  const auto& neutrino = neutrinos->at(bestNuIdx);
283  const LorentzVector w1LVec = lepton.p4()+neutrino.p4();
284 
285  // Continue to hadronic W
286  // We also collect b jets to be used in the next step.
287  dm = 1e9; // Reset dM to very large value.
288  int bestJ1Idx = -1, bestJ2Idx = -1;
289  for ( auto j1Idx = ljetIdxs.begin(); j1Idx != ljetIdxs.end(); ++j1Idx )
290  {
291  const auto& j1 = jets->at(*j1Idx);
292  for ( auto j2Idx = j1Idx+1; j2Idx != ljetIdxs.end(); ++j2Idx )
293  {
294  const auto& j2 = jets->at(*j2Idx);
295  const double mW = (j1.p4()+j2.p4()).mass();
296  if ( mW > 300 ) continue;
297 
298  const double dmNew = abs(mW-wMass_);
299  if ( dmNew < dm )
300  {
301  dm = dmNew;
302  bestJ1Idx = *j1Idx;
303  bestJ2Idx = *j2Idx;
304  }
305  }
306  }
307  if ( bestJ1Idx < 0 ) break;
308  const auto& wJet1 = jets->at(bestJ1Idx);
309  const auto& wJet2 = jets->at(bestJ2Idx);
310  const LorentzVector w2LVec = wJet1.p4() + wJet2.p4();
311 
312  // Now we have leptonic W and hadronic W.
313  // Contiue to top quarks
314  dm = 1e9; // Reset once again for top combination.
315  int bestB1Idx = -1, bestB2Idx = -1;
316  if ( bjetIdxs.size() < 2 ) break;
317  for ( auto b1Idx : bjetIdxs )
318  {
319  const double t1Mass = (w1LVec + jets->at(b1Idx).p4()).mass();
320  if ( t1Mass > 300 ) continue;
321  for ( auto b2Idx : bjetIdxs )
322  {
323  if ( b1Idx == b2Idx ) continue;
324  const double t2Mass = (w2LVec + jets->at(b2Idx).p4()).mass();
325  if ( t2Mass > 300 ) continue;
326 
327  const double dmNew = abs(t1Mass-tMass_) + abs(t2Mass-tMass_);
328  if ( dmNew < dm )
329  {
330  dm = dmNew;
331  bestB1Idx = b1Idx;
332  bestB2Idx = b2Idx;
333  }
334  }
335  }
336  if ( bestB1Idx < 0 ) break;
337  const auto& bJet1 = jets->at(bestB1Idx);
338  const auto& bJet2 = jets->at(bestB2Idx);
339  const LorentzVector t1LVec = w1LVec + bJet1.p4();
340  const LorentzVector t2LVec = w2LVec + bJet2.p4();
341 
342  // Put all of them into candidate collection
343  if ( true ) // Trick to restrict variables' scope to avoid collision
344  {
345  const int lepQ = lepton.charge();
346 
347  reco::GenParticle t1(lepQ*2/3, t1LVec, genVertex_, lepQ*6, 3, false);
348  reco::GenParticle w1(lepQ, w1LVec, genVertex_, lepQ*24, 3, true);
349  reco::GenParticle b1(0, bJet1.p4(), genVertex_, lepQ*5, 1, true);
350  reco::GenParticle l1(lepQ, lepton.p4(), genVertex_, lepton.pdgId(), 1, true);
351  reco::GenParticle n1(0, neutrino.p4(), genVertex_, neutrino.pdgId(), 1, true);
352 
353  reco::GenParticle t2(-lepQ*2/3, t2LVec, genVertex_, -lepQ*6, 3, false);
354  reco::GenParticle w2(-lepQ, w2LVec, genVertex_, -lepQ*24, 3, true);
355  reco::GenParticle b2(0, bJet2.p4(), genVertex_, -lepQ*5, 1, true);
356  reco::GenParticle qA(0, wJet1.p4(), genVertex_, 1, 1, true);
357  reco::GenParticle qB(0, wJet2.p4(), genVertex_, 2, 1, true);
358 
359  pseudoTop->push_back(t1);
360  pseudoTop->push_back(w1);
361  pseudoTop->push_back(b1);
362  pseudoTop->push_back(l1);
363  pseudoTop->push_back(n1);
364 
365  pseudoTop->push_back(t2);
366  pseudoTop->push_back(w2);
367  pseudoTop->push_back(b2);
368  pseudoTop->push_back(qA);
369  pseudoTop->push_back(qB);
370  }
371  }
372  else if ( nLepton == 2 and wCandIdxs.size() >= 2 ) // Dilepton channel
373  {
374  double dm = 1e9;
375  int bestLep1Idx = -1, bestLep2Idx = -1, bestNu1Idx = -1, bestNu2Idx = -1;
376  for ( auto i : wCandIdxs )
377  {
378  const auto& lepton1 = leptons->at(i.first);
379  if ( lepton1.charge() < 0 ) continue;
380  const auto& neutrino1 = neutrinos->at(i.second);
381  const double mW1 = (lepton1.p4()+neutrino1.p4()).mass();
382  if ( mW1 > 300 ) continue;
383  for ( auto j : wCandIdxs )
384  {
385  if ( i == j ) continue;
386  const auto& lepton2 = leptons->at(j.first);
387  if ( lepton2.charge() > 0 ) continue;
388  const auto& neutrino2 = neutrinos->at(j.second);
389  const double mW2 = (lepton2.p4()+neutrino2.p4()).mass();
390  if ( mW2 > 300 ) continue;
391 
392  const double dmNew = abs(mW1-wMass_) + abs(mW2-wMass_);
393  if ( dmNew < dm )
394  {
395  dm = dmNew;
396  bestLep1Idx = i.first;
397  bestLep2Idx = j.first;
398  bestNu1Idx = i.second;
399  bestNu2Idx = j.second;
400  }
401  }
402  }
403  if ( bestLep1Idx < 0 ) break;
404  const auto& lepton1 = leptons->at(bestLep1Idx);
405  const auto& lepton2 = leptons->at(bestLep2Idx);
406  const auto& neutrino1 = neutrinos->at(bestNu1Idx);
407  const auto& neutrino2 = neutrinos->at(bestNu2Idx);
408  const LorentzVector w1LVec = lepton1.p4()+neutrino1.p4();
409  const LorentzVector w2LVec = lepton2.p4()+neutrino2.p4();
410 
411  // Contiue to top quarks
412  dm = 1e9; // Reset once again for top combination.
413  int bestB1Idx = -1, bestB2Idx = -1;
414  if ( bjetIdxs.size() < 2 ) break;
415  for ( auto b1Idx : bjetIdxs )
416  {
417  const double t1Mass = (w1LVec + jets->at(b1Idx).p4()).mass();
418  if ( t1Mass > 300 ) continue;
419  for ( auto b2Idx : bjetIdxs )
420  {
421  if ( b1Idx == b2Idx ) continue;
422  const double t2Mass = (w2LVec + jets->at(b2Idx).p4()).mass();
423  if ( t2Mass > 300 ) continue;
424 
425  const double dmNew = abs(t1Mass-tMass_) + abs(t2Mass-tMass_);
426  if ( dmNew < dm )
427  {
428  dm = dmNew;
429  bestB1Idx = b1Idx;
430  bestB2Idx = b2Idx;
431  }
432  }
433  }
434  if ( bestB1Idx < 0 ) break;
435  const auto& bJet1 = jets->at(bestB1Idx);
436  const auto& bJet2 = jets->at(bestB2Idx);
437  const LorentzVector t1LVec = w1LVec + bJet1.p4();
438  const LorentzVector t2LVec = w2LVec + bJet2.p4();
439 
440  // Put all of them into candidate collection
441  if ( true ) // Trick to restrict variables' scope to avoid collision
442  {
443  const int lep1Q = lepton1.charge();
444  const int lep2Q = lepton2.charge();
445 
446  reco::GenParticle t1(lep1Q*2/3, t1LVec, genVertex_, lep1Q*6, 3, true);
447  reco::GenParticle w1(lep1Q, w1LVec, genVertex_, lep1Q*24, 3, true);
448  reco::GenParticle b1(0, bJet1.p4(), genVertex_, lep1Q*5, 1, true);
449  reco::GenParticle l1(lep1Q, lepton1.p4(), genVertex_, lepton1.pdgId(), 1, true);
450  reco::GenParticle n1(0, neutrino1.p4(), genVertex_, neutrino1.pdgId(), 1, true);
451 
452  reco::GenParticle t2(lep2Q*2/3, t2LVec, genVertex_, lep2Q*6, 3, true);
453  reco::GenParticle w2(lep2Q, w2LVec, genVertex_, lep2Q*24, 3, true);
454  reco::GenParticle b2(0, bJet2.p4(), genVertex_, lep2Q*5, 1, true);
455  reco::GenParticle l2(0, lepton2.p4(), genVertex_, lepton2.pdgId(), 1, true);
456  reco::GenParticle n2(0, neutrino2.p4(), genVertex_, neutrino2.pdgId(), 1, true);
457 
458  pseudoTop->push_back(t1);
459  pseudoTop->push_back(w1);
460  pseudoTop->push_back(b1);
461  pseudoTop->push_back(l1);
462  pseudoTop->push_back(n1);
463 
464  pseudoTop->push_back(t2);
465  pseudoTop->push_back(w2);
466  pseudoTop->push_back(b2);
467  pseudoTop->push_back(l2);
468  pseudoTop->push_back(n2);
469  }
470  }
471 
472  if ( pseudoTop->size() == 10 ) // If pseudtop decay tree is completed
473  {
474  // t->W+b
475  pseudoTop->at(0).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 1)); // t->W
476  pseudoTop->at(0).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 2)); // t->b
477  pseudoTop->at(1).addMother(reco::GenParticleRef(pseudoTopRefHandle, 0)); // t->W
478  pseudoTop->at(2).addMother(reco::GenParticleRef(pseudoTopRefHandle, 0)); // t->b
479 
480  // W->lv or W->jj
481  pseudoTop->at(1).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 3));
482  pseudoTop->at(1).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 4));
483  pseudoTop->at(3).addMother(reco::GenParticleRef(pseudoTopRefHandle, 1));
484  pseudoTop->at(4).addMother(reco::GenParticleRef(pseudoTopRefHandle, 1));
485 
486  // tbar->W-b
487  pseudoTop->at(5).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 6));
488  pseudoTop->at(5).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 7));
489  pseudoTop->at(6).addMother(reco::GenParticleRef(pseudoTopRefHandle, 5));
490  pseudoTop->at(7).addMother(reco::GenParticleRef(pseudoTopRefHandle, 5));
491 
492  // W->jj
493  pseudoTop->at(6).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 8));
494  pseudoTop->at(6).addDaughter(reco::GenParticleRef(pseudoTopRefHandle, 9));
495  pseudoTop->at(8).addMother(reco::GenParticleRef(pseudoTopRefHandle, 6));
496  pseudoTop->at(9).addMother(reco::GenParticleRef(pseudoTopRefHandle, 6));
497  }
498  } while (false);
499 
500  event.put(neutrinos, "neutrinos");
501  event.put(leptons, "leptons");
502  event.put(jets, "jets");
503 
504  event.put(pseudoTop);
505 }
bool isBHadron(const reco::Candidate *p) const
const double leptonMaxEta_
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int i
Definition: DBlmapReader.cc:9
reco::Particle::Point genVertex_
virtual double energy() const =0
energy
const double leptonMinPt_
const double jetMaxEta_
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
edm::EDGetTokenT< edm::View< reco::Candidate > > finalStateToken_
tuple bJet1
do the razor with one or two b jets (medium CSV)
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
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::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
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 double pz() const =0
z coordinate of momentum vector
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
bool isFromHadron(const reco::Candidate *p) const
math::XYZTLorentzVector LorentzVector
std::shared_ptr< JetDef > fjLepDef_
virtual double p() const =0
magnitude of momentum vector
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
bool isnan(float x)
Definition: math.h:13
vector< PseudoJet > jets
tuple dm
Definition: symbols.py:65
bool isNull() const
Checks for null.
Definition: Ptr.h:165
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
Jets made from MC generator particles.
Definition: GenJet.h:24
virtual const Point & vertex() const =0
vertex position
virtual int charge() const =0
electric charge
virtual double py() const =0
y coordinate of momentum vector
edm::EDGetTokenT< edm::View< reco::Candidate > > genParticleToken_
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:169
virtual double px() const =0
x coordinate of momentum vector
virtual int pdgId() const =0
PDG identifier.
const double jetMinPt_
tuple status
Definition: mps_update.py:57
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, edm::EventSetup const &c)
Definition: JetSpecific.cc:41
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector

Member Data Documentation

edm::EDGetTokenT<edm::View<reco::Candidate> > PseudoTopProducer::finalStateToken_
private

Definition at line 31 of file PseudoTopProducer.h.

Referenced by produce(), and PseudoTopProducer().

std::shared_ptr<JetDef> PseudoTopProducer::fjJetDef_
private

Definition at line 37 of file PseudoTopProducer.h.

Referenced by produce(), and PseudoTopProducer().

std::shared_ptr<JetDef> PseudoTopProducer::fjLepDef_
private

Definition at line 37 of file PseudoTopProducer.h.

Referenced by produce(), and PseudoTopProducer().

edm::EDGetTokenT<edm::View<reco::Candidate> > PseudoTopProducer::genParticleToken_
private

Definition at line 32 of file PseudoTopProducer.h.

Referenced by produce(), and PseudoTopProducer().

reco::Particle::Point PseudoTopProducer::genVertex_
private

Definition at line 38 of file PseudoTopProducer.h.

Referenced by produce(), and PseudoTopProducer().

const double PseudoTopProducer::jetMaxEta_
private

Definition at line 33 of file PseudoTopProducer.h.

Referenced by produce().

const double PseudoTopProducer::jetMinPt_
private

Definition at line 33 of file PseudoTopProducer.h.

Referenced by produce().

const double PseudoTopProducer::leptonMaxEta_
private

Definition at line 33 of file PseudoTopProducer.h.

Referenced by produce().

const double PseudoTopProducer::leptonMinPt_
private

Definition at line 33 of file PseudoTopProducer.h.

Referenced by produce().

const double PseudoTopProducer::tMass_
private

Definition at line 34 of file PseudoTopProducer.h.

Referenced by produce().

const double PseudoTopProducer::wMass_
private

Definition at line 34 of file PseudoTopProducer.h.

Referenced by produce().