7 #include <Pythia8/Pythia.h> 8 #include "Pythia8Plugins/HepMC2.h" 14 pythia_->settings.flag(
"ProcessLevel:all",
false);
15 pythia_->settings.flag(
"PartonLevel:FSRinResonances",
false);
16 pythia_->settings.flag(
"ProcessLevel:resonanceDecays",
false);
23 while (pdt.nextId(pid) > pid) {
24 pid = pdt.nextId(pid);
25 pdt.mayDecay(pid,
false);
30 std::vector<std::unique_ptr<fastsim::Particle> >& secondaries,
31 CLHEP::HepRandomEngine& engine)
const {
36 int pid = particle.
pdgId();
44 pythia_->event.reset();
47 Pythia8::Particle pythiaParticle(pid,
62 pythia_->event.append(pythiaParticle);
64 int nentries_before = pythia_->event.size();
66 pythia_->particleData.mayDecay(pid,
true);
70 pythia_->particleData.mayDecay(pid,
false);
71 int nentries_after = pythia_->event.size();
73 if (nentries_after <= nentries_before)
77 for (
int ipart = nentries_before; ipart < nentries_after; ipart++) {
78 Pythia8::Particle& daughter = pythia_->event[ipart];
87 secondaries.back()->setMotherDeltaR(particle.
momentum());
90 secondaries.back()->setMotherSimTrackIndex(particle.
simTrackIndex());
std::unique_ptr< gen::P8RndmEngine > pythiaRandomEngine_
Instance of pythia Random Engine.
int getMotherPdgId() const
Get pdgIdto mother particle (only makes sense if mother and daughter charged).
void decay(const Particle &particle, std::vector< std::unique_ptr< Particle > > &secondaries, CLHEP::HepRandomEngine &engine) const
Decay particle using pythia.
const math::XYZTLorentzVector & position() const
Return position of the particle.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
~Decayer()
Default destructor.
bool isExotic(int pdgid_)
Abs< T >::type abs(const T &t)
int simTrackIndex() const
Return index of the SimTrack.
HepPDT::ParticleData ParticleData
int pdgId() const
Return pdgId of the particle.
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
double charge() const
Return charge of the particle.
Decayer()
Default Constructor.
double getMotherDeltaR() const
Get delta R to mother particle (only makes sense if mother and daughter charged). ...
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
std::unique_ptr< Pythia8::Pythia > pythia_
Instance of pythia.