12#include <TClonesArray.h>
28#ifdef __GENIE_PYTHIA6_ENABLED__
29#if ROOT_VERSION_CODE >= ROOT_VERSION(5,15,6)
30#include <TMCParticle.h>
32#include <TMCParticle6.h>
40extern "C" void py2ent_(
int *,
int *,
int *,
double *);
61 #ifdef __GENIE_PYTHIA6_ENABLED__
66#ifdef __GENIE_PYTHIA6_ENABLED__
70 <<
"Calling GENIE/PYTHIA6 hadronization modules without enabling PYTHIA6";
77#ifdef __GENIE_PYTHIA6_ENABLED__
82#ifdef __GENIE_PYTHIA6_ENABLED__
84 LOG(
"Pythia6Had",
pNOTICE) <<
"Running PYTHIA6 hadronizer";
88 double W = kinematics.
W();
100 fPythia->GetPrimaries();
101 TClonesArray * pythia_particles =
102 (TClonesArray *) fPythia->ImportParticles(
"All");
107 int np = pythia_particles->GetEntries();
111 TLorentzVector p4Had = kinematics.
HadSystP4();
114 TVector3 unitvq = p4Had.Vect().Unit();
117 TVector3 beta(0,0,p4Had.P()/p4Had.Energy());
124 int mom =
event->FinalStateHadronicSystemPosition();
129 const TLorentzVector & vtx = *(neutrino->
X4());
134 TIter particle_iter(pythia_particles);
135 while( (p = (TMCParticle *) particle_iter.Next()) ) {
137 int particle_pdg_code = p->GetKF();
138 int pythia_particle_status = p->GetKS();
141 if(pythia_particle_status == 1) {
146 <<
"Hadronization failed! Bare quarks appear in final state!";
156 TLorentzVector p4o(p->GetPx(), p->GetPy(), p->GetPz(), p->GetEnergy());
158 TVector3 p3 = p4o.Vect();
160 TLorentzVector p4(p3,p4o.Energy());
168 bool is_gamma = (particle_pdg_code ==
kPdgGamma);
171 bool not_hadr = is_gamma || is_nu || is_lchg;
175 int mother1 = mom + p->GetParent();
177 int daughter1 = (p->GetFirstChild() <= 0 ) ? -1 : mom + p->GetFirstChild();
178 int daughter2 = (p->GetLastChild() <= 0 ) ? -1 : mom + p->GetLastChild();
199 <<
"Adding final state particle pdgc = " << particle.
Pdg()
200 <<
" with status = " << particle.
Status();
203 event->AddParticle(particle);
214#ifdef __GENIE_PYTHIA6_ENABLED__
225#ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
227 <<
"Original PYTHIA6 decay flags:"
244#ifdef __GENIE_PYTHIA6_ENABLED__
246 fPythia->SetMDCY(fPythia->Pycomp(
kPdgPi0),
248 fPythia->SetMDCY(fPythia->Pycomp(
kPdgK0),
270#ifdef __GENIE_PYTHIA6_ENABLED__
272fPythia->SetMDCY(fPythia->Pycomp(
kPdgPi0),
274fPythia->SetMDCY(fPythia->Pycomp(
kPdgK0),
310#ifdef __GENIE_PYTHIA6_ENABLED__
318 fPythia->SetPARJ(41,
fLunda );
319 fPythia->SetPARJ(42,
fLundb );
329#ifdef __GENIE_PYTHIA6_ENABLED__
330 fPythia = TPythia6::Instance();
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
void py2ent_(int *, int *, int *, double *)
virtual const Registry & GetConfig(void) const
virtual void Configure(const Registry &config)
STDHEP-like event record entry that can fit a particle or a nucleus.
const TLorentzVector * X4(void) const
GHepStatus_t Status(void) const
GENIE's GHEP MC event record.
const Target & Tgt(void) const
Summary information for an interaction.
const Kinematics & Kine(void) const
const InitialState & InitState(void) const
Generated/set kinematical variables for an event.
const TLorentzVector & HadSystP4(void) const
double W(bool selected=false) const
virtual ~Pythia6Hadro2019()
void Configure(const Registry &config)
void SetDesiredDecayFlags(void) const
void RestoreOriginalDecayFlags(void) const
bool Hadronize(GHepRecord *event) const
void ProcessEventRecord(GHepRecord *event) const
void CopyOriginalDecayFlags(void) const
double fLundb
Lund b parameter.
double fSSBarSuppression
ssbar suppression
double fNonGaussianPt2Tail
non gaussian pt2 tail parameterization
double fDiQuarkSuppression
di-quark suppression parameter
double fLundaDiq
adjustment of Lund a for di-quark
double fGaussianPt2
gaussian pt2 distribution width
virtual void LoadConfig(void)
virtual void Initialize(void)
double fRemainingECutoff
remaining E cutoff stopping fragmentation
double fLunda
Lund a parameter.
double fSVMesonSuppression
strange vector meson suppression
virtual void ProcessEventRecord(GHepRecord *event) const
double fLightVMesonSuppression
light vector meson suppression
static RandomGen * Instance()
Access instance.
A registry. Provides the container for algorithm configuration parameters.
bool IsNucleus(void) const
bool IsChargedLepton(int pdgc)
bool IsNeutralLepton(int pdgc)
THE MAIN GENIE PROJECT NAMESPACE
@ kIStDISPreFragmHadronicState
const int kPdgP33m1232_Delta0
const int kPdgP33m1232_DeltaPP
enum genie::EGHepStatus GHepStatus_t
const int kPdgP33m1232_DeltaM
const int kPdgP33m1232_DeltaP