CMS 3D CMS Logo

List of all members | Public Member Functions
Pythia8PartonSelector Class Reference

Pythia8 parton selector derived from the base parton selector. More...

#include <PhysicsTools/JetMCAlgos/interface/Pythia8PartonSelector.h>

Inheritance diagram for Pythia8PartonSelector:
BasePartonSelector

Public Member Functions

 Pythia8PartonSelector ()
 
void run (const edm::Handle< reco::GenParticleCollection > &particles, std::unique_ptr< reco::GenParticleRefVector > &partons) override
 
 ~Pythia8PartonSelector () override
 
- Public Member Functions inherited from BasePartonSelector
 BasePartonSelector ()
 
virtual ~BasePartonSelector ()
 

Detailed Description

Pythia8 parton selector derived from the base parton selector.

Definition at line 10 of file Pythia8PartonSelector.h.

Constructor & Destructor Documentation

◆ Pythia8PartonSelector()

Pythia8PartonSelector::Pythia8PartonSelector ( )

This is a Pythia8-specific parton selector that selects all partons that don't have other partons as daughters, i.e., partons from the end of the parton showering sequence. An explanation of the particle status codes returned by Pythia8 can be found in Pythia8 online manual (http://home.thep.lu.se/~torbjorn/pythia81html/ParticleProperties.html).

Definition at line 11 of file Pythia8PartonSelector.cc.

11 {}

◆ ~Pythia8PartonSelector()

Pythia8PartonSelector::~Pythia8PartonSelector ( )
override

Definition at line 13 of file Pythia8PartonSelector.cc.

13 {}

Member Function Documentation

◆ run()

void Pythia8PartonSelector::run ( const edm::Handle< reco::GenParticleCollection > &  particles,
std::unique_ptr< reco::GenParticleRefVector > &  partons 
)
overridevirtual

Reimplemented from BasePartonSelector.

Definition at line 15 of file Pythia8PartonSelector.cc.

16  {
17  // loop over particles and select partons
18  for (reco::GenParticleCollection::const_iterator it = particles->begin(); it != particles->end(); ++it) {
19  int status = it->status();
20  if (status == 1)
21  continue; // skip stable particles
22  if (status == 2)
23  continue; // skip decayed Standard Model hadrons and leptons
24  if (!CandMCTagUtils::isParton(*it))
25  continue; // skip particle if not a parton
26 
27  // check if the parton has other partons as daughters
28  int nparton_daughters = 0;
29  for (unsigned i = 0; i < it->numberOfDaughters(); ++i) {
30  if (CandMCTagUtils::isParton(*(it->daughter(i))))
31  ++nparton_daughters;
32  }
33 
34  if (nparton_daughters == 0)
35  partons->push_back(reco::GenParticleRef(particles, it - particles->begin()));
36  }
37 
38  return;
39 }

References mps_fire::i, CandMCTagUtils::isParton(), ecalTrigSettings_cff::particles, dqmAnalyzer_cff::partons, and mps_update::status.

mps_fire.i
i
Definition: mps_fire.py:428
mps_update.status
status
Definition: mps_update.py:68
CandMCTagUtils::isParton
bool isParton(const reco::Candidate &c)
Definition: CandMCTag.cc:46
ecalTrigSettings_cff.particles
particles
Definition: ecalTrigSettings_cff.py:11
edm::Ref< GenParticleCollection >
dqmAnalyzer_cff.partons
partons
Definition: dqmAnalyzer_cff.py:28