#include <PythiaDecays.h>
Public Member Functions | |
const DaughterParticleList & | particleDaughtersPy6 (ParticlePropagator &particle) |
const DaughterParticleList & | particleDaughtersPy8 (ParticlePropagator &particle) |
PythiaDecays (std::string program) | |
~PythiaDecays () | |
Private Attributes | |
std::auto_ptr< Pythia8::Pythia > | decayer |
std::string | program_ |
Pythia6jets * | pyjets |
gen::Pythia6Service * | pyservice |
std::auto_ptr< Pythia8::Pythia > | pythia |
DaughterParticleList | theList |
Definition at line 27 of file PythiaDecays.h.
PythiaDecays::PythiaDecays | ( | std::string | program | ) |
Definition at line 21 of file PythiaDecays.cc.
References gather_cfg::cout, decayer, program_, pyjets, pyservice, pythia, and Pythia6jets::Pythia6jets().
{ program_=program; if (program_ == "pythia6") { pyjets = new Pythia6jets(); pyservice = new gen::Pythia6Service(); // The PYTHIA decay tables will be initialized later } else if (program_ == "pythia8") { pythia.reset(new Pythia8::Pythia); decayer.reset(new Pythia8::Pythia); RandomP8* RP8 = new RandomP8(); pythia->setRndmEnginePtr(RP8); decayer->setRndmEnginePtr(RP8); } else { std::cout << "WARNING: you are requesting an option which is not available in PythiaDecays::PythiaDecays " << std::endl; } }
PythiaDecays::~PythiaDecays | ( | ) |
const DaughterParticleList & PythiaDecays::particleDaughtersPy6 | ( | ParticlePropagator & | particle | ) |
Definition at line 108 of file PythiaDecays.cc.
References i, Pythia6jets::k(), RawParticle::mass(), max(), Pythia6jets::n(), Pythia6jets::p(), RawParticle::pid(), pyjets, pyservice, PYTHIA6PYDECY, RawParticle::T(), theList, Pythia6jets::v(), RawParticle::X(), RawParticle::Y(), and RawParticle::Z().
Referenced by TrajectoryManager::updateWithDaughters().
{ gen::Pythia6Service::InstanceWrapper guard(pyservice); // grab Py6 context // Pythia6jets pyjets; int ip; pyjets->k(1,1) = 1; pyjets->k(1,2) = particle.pid(); pyjets->p(1,1) = particle.Px(); pyjets->p(1,2) = particle.Py(); pyjets->p(1,3) = particle.Pz(); pyjets->p(1,4) = std::max(particle.mass(),particle.e()); pyjets->p(1,5) = particle.mass(); pyjets->v(1,1) = particle.X(); pyjets->v(1,2) = particle.Y(); pyjets->v(1,3) = particle.Z(); pyjets->v(1,4) = particle.T(); pyjets->n() = 1; ip = 1; PYTHIA6PYDECY(&ip); // Fill the list of daughters theList.clear(); if ( pyjets->n()==1 ) return theList; theList.resize(pyjets->n()-1,RawParticle()); for (int i=2;i<=pyjets->n();++i) { theList[i-2].SetXYZT(pyjets->p(i,1),pyjets->p(i,2), pyjets->p(i,3),pyjets->p(i,4)); theList[i-2].setVertex(pyjets->v(i,1),pyjets->v(i,2), pyjets->v(i,3),pyjets->v(i,4)); theList[i-2].setID(pyjets->k(i,2)); theList[i-2].setMass(pyjets->p(i,5)); } return theList; }
const DaughterParticleList & PythiaDecays::particleDaughtersPy8 | ( | ParticlePropagator & | particle | ) |
Definition at line 51 of file PythiaDecays.cc.
References decayer, RawParticle::mass(), RawParticle::momentum(), RawParticle::pid(), RawParticle::T(), theList, RawParticle::X(), RawParticle::Y(), and RawParticle::Z().
Referenced by TrajectoryManager::updateWithDaughters().
{ theList.clear(); /* std::cout << "###" << std::endl; std::cout << "particle.PDGname() == " << particle.PDGname() << std::endl; std::cout << "particle.status() == " << particle.status() << std::endl; std::cout << "particle.pid() == " << particle.pid() << std::endl; std::cout << "particle.momentum().x() == " << particle.momentum().x() << std::endl; std::cout << "particle.Px() == " << particle.Px() << std::endl; std::cout << "particle.X() == " << particle.X() << std::endl; std::cout << "particle.mass() == " << particle.mass() << std::endl; std::cout << "particle.PDGmass() == " << particle.PDGmass() << std::endl; std::cout << "particle.PDGcTau() == " << particle.PDGcTau() << std::endl; std::cout << "(decayer->particleData).tau0( particle.pid() ) == " << (decayer->particleData).tau0( particle.pid() ) << std::endl; */ // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Pythia8Hadronizer.cc decayer->event.reset(); Pythia8::Particle py8part( particle.pid(), 93, 0, 0, 0, 0, 0, 0, particle.momentum().x(), // note: momentum().x() and Px() are the same particle.momentum().y(), particle.momentum().z(), particle.momentum().t(), particle.mass() ); py8part.vProd( particle.X(), particle.Y(), particle.Z(), particle.T() ); py8part.tau( (decayer->particleData).tau0( particle.pid() ) ); decayer->event.append( py8part ); int nentries = decayer->event.size(); if ( !decayer->event[nentries-1].mayDecay() ) return theList; decayer->next(); int nentries1 = decayer->event.size(); if ( nentries1 <= nentries ) return theList; //same number of particles, no decays... Pythia8::Particle& py8daughter = decayer->event[nentries]; // the 1st daughter // DO I NEED THIS LINE? theList.resize(nentries,RawParticle()); for ( int ipart=nentries+1; ipart<nentries1; ipart++ ) { py8daughter = decayer->event[ipart]; theList[ipart-nentries-1].SetXYZT( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() ); theList[ipart-nentries-1].setVertex( py8daughter.xProd(), py8daughter.yProd(), py8daughter.zProd(), py8daughter.tProd() ); theList[ipart-nentries-1].setID( py8daughter.id() ); theList[ipart-nentries-1].setMass( py8daughter.m() ); } return theList; }
std::auto_ptr<Pythia8::Pythia> PythiaDecays::decayer [private] |
Definition at line 46 of file PythiaDecays.h.
Referenced by particleDaughtersPy8(), and PythiaDecays().
std::string PythiaDecays::program_ [private] |
Definition at line 40 of file PythiaDecays.h.
Referenced by PythiaDecays(), and ~PythiaDecays().
Pythia6jets* PythiaDecays::pyjets [private] |
Definition at line 43 of file PythiaDecays.h.
Referenced by particleDaughtersPy6(), PythiaDecays(), and ~PythiaDecays().
gen::Pythia6Service* PythiaDecays::pyservice [private] |
Definition at line 42 of file PythiaDecays.h.
Referenced by particleDaughtersPy6(), PythiaDecays(), and ~PythiaDecays().
std::auto_ptr<Pythia8::Pythia> PythiaDecays::pythia [private] |
Definition at line 45 of file PythiaDecays.h.
Referenced by PythiaDecays().
DaughterParticleList PythiaDecays::theList [private] |
Definition at line 39 of file PythiaDecays.h.
Referenced by particleDaughtersPy6(), and particleDaughtersPy8().