GENIEGenerator
Loading...
Searching...
No Matches
LeptoHadPythia8.h
Go to the documentation of this file.
1//____________________________________________________________________________
2/*!
3
4\class genie::LeptoHadPythia8
5
6\brief Provides access to the PYTHIA8 hadronization. \n
7
8\author Alfonso Garcia <aagarciasoto \at km3net.de>
9 IFIC (Valencia)
10
11\created December 12, 2024
12
13\cpright Copyright (c) 2003-2025, The GENIE Collaboration
14 For the full text of the license visit http://copyright.genie-mc.org
15 or see $GENIE/LICENSE
16*/
17//____________________________________________________________________________
18
19#ifndef _LEPTO_HAD_PYTHIA8__H_
20#define _LEPTO_HAD_PYTHIA8__H_
21
33
34#ifdef __GENIE_PYTHIA8_ENABLED__
36#endif
37
38using namespace genie;
39using namespace genie::constants;
40using namespace genie::utils::math;
41
42namespace genie {
43
44class GHepParticle;
45
47
48public:
50 LeptoHadPythia8(string config);
51 virtual ~LeptoHadPythia8();
52
53 //-- implement the HadronizationModelI interface
54 void ProcessEventRecord(GHepRecord * event) const;
55
56 //-- overload the Algorithm::Configure() methods to load private data
57 // members from configuration options
58 void Configure(const Registry & config);
59 void Configure(string config);
60
61private:
62
63 bool Hadronize (GHepRecord * event) const;
64 void Initialize (void) const;
65 void LoadConfig (void);
66
67 int getMeson (int,int,double) const; ///< create meson
68 int getBaryon (int,int,double) const; ///< create baryon
69 double getRandomZ (double,double) const; ///< fragmentation
70 double Afrag; ///< fragmentation parameter
71 double Bfrag; ///< fragmentation parameter
72 double mesonRateSum[3]; ///<meson parameter
73 double CGOct[6] = { 0.75, 0.5, 0., 0.1667, 0.0833, 0.1667}; ///<baryon parameter
74 double CGDec[6] = { 0.00, 0.0, 1., 0.3333, 0.6667, 0.3333}; ///<baryon parameter
75 double CGSum[6]; ///< baryon parameter
76
77 //-- configuration parameters
78 int fMaxIterHad; // Maxmium number of iterations to look for a combination of hadrons
79 double fPrimordialKT; // Width of Gaussian distribution for the primordial transverse momentum kT of partons in the nucleon.
80 double fRemnantPT; // Width of Gaussian distribution in transverse momentum when a non-trivial target remnant is split into two particles
81 double fMinESinglet; // It is, with quark masses added, used to define the minimum allowable energy of a colour-singlet parton system.
82 double fWmin; // Minimum value of W
83
84 // PYTHIA physics configuration parameters used
85 double fSSBarSuppression; ///< ssbar suppression
86 double fGaussianPt2; ///< gaussian pt2 distribution width
87 double fNonGaussianPt2Tail; ///< non gaussian pt2 tail parameterization
88 double fRemainingECutoff; ///< remaining E cutoff stopping fragmentation
89 double fDiQuarkSuppression; ///< di-quark suppression parameter
90 double fLightVMesonSuppression; ///< light vector meson suppression
91 double fSVMesonSuppression; ///< strange vector meson suppression
92 double fLunda; ///< Lund a parameter
93 double fLundb; ///< Lund b parameter
94 double fLundaDiq; ///< adjustment of Lund a for di-quark
95
96#ifdef __GENIE_PYTHIA8_ENABLED__
97 mutable Pythia8::Pythia * fPythia; ///< PYTHIA8 object
98#endif
99
100};
101
102} // genie namespace
103
104#endif // _LEPTO_HAD_PYTHIA8__H_
105
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
STDHEP-like event record entry that can fit a particle or a nucleus.
GENIE's GHEP MC event record.
Definition GHepRecord.h:45
double mesonRateSum[3]
meson parameter
void Configure(const Registry &config)
double fSSBarSuppression
ssbar suppression
double fNonGaussianPt2Tail
non gaussian pt2 tail parameterization
double fLundaDiq
adjustment of Lund a for di-quark
double fLunda
Lund a parameter.
int getBaryon(int, int, double) const
create baryon
double getRandomZ(double, double) const
fragmentation
double CGOct[6]
baryon parameter
double CGDec[6]
baryon parameter
void ProcessEventRecord(GHepRecord *event) const
int getMeson(int, int, double) const
create meson
double fGaussianPt2
gaussian pt2 distribution width
double Bfrag
fragmentation parameter
double fDiQuarkSuppression
di-quark suppression parameter
double fLightVMesonSuppression
light vector meson suppression
double fRemainingECutoff
remaining E cutoff stopping fragmentation
void Initialize(void) const
double CGSum[6]
baryon parameter
bool Hadronize(GHepRecord *event) const
double fSVMesonSuppression
strange vector meson suppression
double Afrag
fragmentation parameter
double fLundb
Lund b parameter.
A registry. Provides the container for algorithm configuration parameters.
Definition Registry.h:65
Basic constants.
Simple functions for loading and reading nucleus dependent keys from config files.
Simple mathematical utilities not found in ROOT's TMath.
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25