GENIEGenerator
Loading...
Searching...
No Matches
genie::LeptoHadPythia6 Class Reference

Provides access to the PYTHIA6 hadronization.
. More...

#include <LeptoHadPythia6.h>

Inheritance diagram for genie::LeptoHadPythia6:
[legend]
Collaboration diagram for genie::LeptoHadPythia6:
[legend]

Public Member Functions

 LeptoHadPythia6 ()
 LeptoHadPythia6 (string config)
virtual ~LeptoHadPythia6 ()
void ProcessEventRecord (GHepRecord *event) const
void Configure (const Registry &config)
void Configure (string config)
Public Member Functions inherited from genie::EventRecordVisitorI
virtual ~EventRecordVisitorI ()
Public Member Functions inherited from genie::Algorithm
virtual ~Algorithm ()
virtual void FindConfig (void)
virtual const RegistryGetConfig (void) const
RegistryGetOwnedConfig (void)
virtual const AlgIdId (void) const
 Get algorithm ID.
virtual AlgStatus_t GetStatus (void) const
 Get algorithm status.
virtual bool AllowReconfig (void) const
virtual AlgCmp_t Compare (const Algorithm *alg) const
 Compare with input algorithm.
virtual void SetId (const AlgId &id)
 Set algorithm ID.
virtual void SetId (string name, string config)
const AlgorithmSubAlg (const RgKey &registry_key) const
void AdoptConfig (void)
void AdoptSubstructure (void)
virtual void Print (ostream &stream) const
 Print algorithm info.

Private Member Functions

bool Hadronize (GHepRecord *event) const
void Initialize (void) const
void LoadConfig (void)

Private Attributes

int fMaxIterHad
double fPrimordialKT
double fRemnantPT
double fMinESinglet
double fWmin
double fSSBarSuppression
 ssbar suppression
double fGaussianPt2
 gaussian pt2 distribution width
double fNonGaussianPt2Tail
 non gaussian pt2 tail parameterization
double fRemainingECutoff
 remaining E cutoff stopping fragmentation
double fDiQuarkSuppression
 di-quark suppression parameter
double fLightVMesonSuppression
 light vector meson suppression
double fSVMesonSuppression
 strange vector meson suppression
double fLunda
 Lund a parameter.
double fLundb
 Lund b parameter.
double fLundaDiq
 adjustment of Lund a for di-quark

Additional Inherited Members

Static Public Member Functions inherited from genie::Algorithm
static string BuildParamVectKey (const std::string &comm_name, unsigned int i)
static string BuildParamVectSizeKey (const std::string &comm_name)
static string BuildParamMatKey (const std::string &comm_name, unsigned int i, unsigned int j)
static string BuildParamMatRowSizeKey (const std::string &comm_name)
static string BuildParamMatColSizeKey (const std::string &comm_name)
Protected Member Functions inherited from genie::EventRecordVisitorI
 EventRecordVisitorI ()
 EventRecordVisitorI (string name)
 EventRecordVisitorI (string name, string config)
Protected Member Functions inherited from genie::Algorithm
 Algorithm ()
 Algorithm (string name)
 Algorithm (string name, string config)
void Initialize (void)
void DeleteConfig (void)
void DeleteSubstructure (void)
RegistryExtractLocalConfig (const Registry &in) const
RegistryExtractLowerConfig (const Registry &in, const string &alg_key) const
 Split an incoming configuration Registry into a block valid for the sub-algo identified by alg_key.
template<class T>
bool GetParam (const RgKey &name, T &p, bool is_top_call=true) const
template<class T>
bool GetParamDef (const RgKey &name, T &p, const T &def) const
template<class T>
int GetParamVect (const std::string &comm_name, std::vector< T > &v, bool is_top_call=true) const
 Handle to load vectors of parameters.
int GetParamVectKeys (const std::string &comm_name, std::vector< RgKey > &k, bool is_top_call=true) const
template<class T>
int GetParamMat (const std::string &comm_name, TMatrixT< T > &mat, bool is_top_call=true) const
 Handle to load matrix of parameters.
template<class T>
int GetParamMatSym (const std::string &comm_name, TMatrixTSym< T > &mat, bool is_top_call=true) const
int GetParamMatKeys (const std::string &comm_name, std::vector< RgKey > &k, bool is_top_call=true) const
int AddTopRegistry (Registry *rp, bool owns=true)
 add registry with top priority, also update ownership
int AddLowRegistry (Registry *rp, bool owns=true)
 add registry with lowest priority, also update ownership
int MergeTopRegistry (const Registry &r)
int AddTopRegisties (const vector< Registry * > &rs, bool owns=false)
 Add registries with top priority, also udated Ownerships.
Protected Attributes inherited from genie::Algorithm
bool fAllowReconfig
bool fOwnsSubstruc
 true if it owns its substructure (sub-algs,...)
AlgId fID
 algorithm name and configuration set
vector< Registry * > fConfVect
vector< bool > fOwnerships
 ownership for every registry in fConfVect
AlgStatus_t fStatus
 algorithm execution status
AlgMapfOwnedSubAlgMp
 local pool for owned sub-algs (taken out of the factory pool)

Detailed Description

Provides access to the PYTHIA6 hadronization.
.

Author
Alfonso Garcia <aagarciasoto \at km3net.de> IFIC (Valencia)
Created:\n December 12, 2024
License:\n Copyright (c) 2003-2025, The GENIE Collaboration
For the full text of the license visit http://copyright.genie-mc.org or see $GENIE/LICENSE

Definition at line 46 of file LeptoHadPythia6.h.

Constructor & Destructor Documentation

◆ LeptoHadPythia6() [1/2]

LeptoHadPythia6::LeptoHadPythia6 ( )

Definition at line 31 of file LeptoHadPythia6.cxx.

31 :
32EventRecordVisitorI("genie::LeptoHadPythia6")
33{
34 this->Initialize();
35}
void Initialize(void) const

References genie::EventRecordVisitorI::EventRecordVisitorI(), and Initialize().

◆ LeptoHadPythia6() [2/2]

LeptoHadPythia6::LeptoHadPythia6 ( string config)

Definition at line 37 of file LeptoHadPythia6.cxx.

37 :
38EventRecordVisitorI("genie::LeptoHadPythia6", config)
39{
40 this->Initialize();
41}

References genie::EventRecordVisitorI::EventRecordVisitorI(), and Initialize().

◆ ~LeptoHadPythia6()

LeptoHadPythia6::~LeptoHadPythia6 ( )
virtual

Definition at line 43 of file LeptoHadPythia6.cxx.

44{
45
46}

Member Function Documentation

◆ Configure() [1/2]

void LeptoHadPythia6::Configure ( const Registry & config)
virtual

Configure the algorithm with an external registry The registry is merged with the top level registry if it is owned, Otherwise a copy of it is added with the highest priority

Reimplemented from genie::Algorithm.

Definition at line 429 of file LeptoHadPythia6.cxx.

430{
431 Algorithm::Configure(config);
432 this->LoadConfig();
433}
virtual void Configure(const Registry &config)
Definition Algorithm.cxx:62

References genie::Algorithm::Configure(), and LoadConfig().

◆ Configure() [2/2]

void LeptoHadPythia6::Configure ( string config)
virtual

Configure the algorithm from the AlgoConfigPool based on param_set string given in input An algorithm contains a vector of registries coming from different xml configuration files, which are loaded according a very precise prioriy This methods will load a number registries in order of priority: 1) "Tunable" parameter set from CommonParametes. This is loaded with the highest prioriry and it is designed to be used for tuning procedure Usage not expected from the user. 2) For every string defined in "CommonParame" the corresponding parameter set will be loaded from CommonParameter.xml 3) parameter set specified by the config string and defined in the xml file of the algorithm 4) if config is not "Default" also the Default parameter set from the same xml file will be loaded Effectively this avoids the repetion of a parameter when it is not changed in the requested configuration

Reimplemented from genie::Algorithm.

Definition at line 435 of file LeptoHadPythia6.cxx.

436{
437 Algorithm::Configure(config);
438 this->LoadConfig();
439}

References genie::Algorithm::Configure(), and LoadConfig().

◆ Hadronize()

bool LeptoHadPythia6::Hadronize ( GHepRecord * event) const
private

Definition at line 64 of file LeptoHadPythia6.cxx.

69{
70
71#ifdef __GENIE_PYTHIA6_ENABLED__
72 // Compute kinematics of hadronic system with energy/momentum conservation
73 LongLorentzVector p4v( * event->Probe()->P4() );
74 LongLorentzVector p4N( * event->HitNucleon()->P4() );
75 LongLorentzVector p4l( * event->FinalStatePrimaryLepton()->P4() );
76 LongLorentzVector p4Hadlong( p4v.Px()+p4N.Px()-p4l.Px(), p4v.Py()+p4N.Py()-p4l.Py(), p4v.Pz()+p4N.Pz()-p4l.Pz(), p4v.E()+p4N.E()-p4l.E() );
77
78 LOG("LeptoHad", pDEBUG) << "v [LAB']: " << p4v.E() << " // " << p4v.M2() << " // [ " << p4v.Dx() << " , " << p4v.Dy() << " , " << p4v.Dz() << " ]";
79 LOG("LeptoHad", pDEBUG) << "N [LAB']: " << p4N.E() << " // " << p4N.M2() << " // [ " << p4N.Dx() << " , " << p4N.Dy() << " , " << p4N.Dz() << " ]";
80 LOG("LeptoHad", pDEBUG) << "l [LAB']: " << p4l.E() << " // " << p4l.M2() << " // [ " << p4l.Dx() << " , " << p4l.Dy() << " , " << p4l.Dz() << " ]";
81 LOG("LeptoHad", pDEBUG) << "H [LAB']: " << p4Hadlong.E() << " // " << p4Hadlong.M2() << " // [ " << p4Hadlong.Dx() << " , " << p4Hadlong.Dy() << " , " << p4Hadlong.Dz() << " ]";
82
83 // Translate from long double to double
84 const TLorentzVector & vtx = *( event->Probe()->X4());
85 TLorentzVector p4Had( (double)p4Hadlong.Px(), (double)p4Hadlong.Py(), (double)p4Hadlong.Pz(), (double)p4Hadlong.E() );
86 event->AddParticle(kPdgHadronicSyst, kIStDISPreFragmHadronicState, event->HitNucleonPosition(),-1,-1,-1, p4Had, vtx);
87
88 const Interaction * interaction = event->Summary();
89 interaction->KinePtr()->SetHadSystP4(p4Had);
90 interaction->KinePtr()->SetW(p4Hadlong.M());
91
92 double W = interaction->Kine().W();
93 if(W < fWmin) {
94 LOG("LeptoHad", pWARN) << "Low invariant mass, W = " << W << " GeV!!";
95 return false;
96 }
97
98 const XclsTag & xclstag = interaction->ExclTag();
99 const Target & target = interaction->InitState().Tgt();
100
101 assert(target.HitQrkIsSet());
102
103 bool isp = pdg::IsProton(target.HitNucPdg());
104 int hit_quark = target.HitQrkPdg();
105 int frag_quark = xclstag.FinalQuarkPdg();
106
107 LOG("LeptoHad", pDEBUG) << "Hit nucleon pdgc = " << target.HitNucPdg() << ", W = " << W;
108 LOG("LeptoHad", pDEBUG) << "Selected hit quark pdgc = " << hit_quark << " // Fragmentation quark = " << frag_quark;
109
110 RandomGen * rnd = RandomGen::Instance();
111
112 //
113 // Generate the hadron combination to input PYTHIA
114 //
115
116 //If the hit quark is a d we have these options:
117 /* uud(->q) => uu + q */
118 /* uud d(->q)db => uu + q (d valence and db sea annihilates)*/
119 /* udd(->q) => ud + q */
120 /* udd d(->q)db => ud + q (d valence and db sea annihilates)*/
121 if ( pdg::IsDQuark(hit_quark) ) {
122 // choose diquark system depending on proton or neutron
123 int diquark = 0;
124 if (isp) diquark = kPdgUUDiquarkS1;
125 else diquark = rnd->RndHadro().Rndm()>0.75 ? kPdgUDDiquarkS1 : kPdgUDDiquarkS0;
126 // Check that the trasnferred energy is higher than the mass of the produced quarks
127 double m_frag = PDGLibrary::Instance()->Find(frag_quark)->Mass();
128 double m_diquark = PDGLibrary::Instance()->Find(diquark)->Mass();
129 if( W <= m_frag + m_diquark + fMinESinglet ) {
130 LOG("LeptoHad", pWARN) << "Low invariant mass, W = " << W << " GeV! Returning a null list";
131 LOG("LeptoHad", pWARN) << "frag_quark = " << frag_quark << " -> m = " << m_frag;
132 LOG("LeptoHad", pWARN) << "diquark = " << diquark << " -> m = " << m_diquark;
133 return 0;
134 }
135 // Input the two particles to PYTHIA back to back in the CM frame
136 double e_frag = (W*W - m_diquark*m_diquark + m_frag*m_frag)/2./W;
137 double e_diquark = (W*W + m_diquark*m_diquark - m_frag*m_frag)/2./W;
138 fPythia->Py1ent( -1, frag_quark, e_frag, 0., 0. ); //k(1,2) = 2
139 // If a top quark is produced we decay it because it does not hadronize
140 if ( pdg::IsTQuark(frag_quark) ) {
141 int ip = 1;
142 pydecy_(&ip);
143 }
144 fPythia->Py1ent( fPythia->GetN()+1, diquark, e_diquark, fPythia->GetPARU(1), 0. ); //k(2,2) = 1
145 }
146
147 //If the hit quark is a u we have these options:
148 /* u(->q)ud => ud + q */
149 /* uud u(->q)ub => ud + q (u valence and ub sea annihilates)*/
150 /* u(->q)dd => dd + q */
151 /* udd u(->q)ub => dd + q (u valence and ub sea annihilates)*/
152 else if ( pdg::IsUQuark(hit_quark) ) {
153 // choose diquark system depending on proton or neutron
154 int diquark = 0;
155 if (isp) diquark = rnd->RndHadro().Rndm()>0.75 ? kPdgUDDiquarkS1 : kPdgUDDiquarkS0;
156 else diquark = kPdgDDDiquarkS1;
157 // Check that the trasnferred energy is higher than the mass of the produced quarks.
158 double m_frag = PDGLibrary::Instance()->Find(frag_quark)->Mass();
159 double m_diquark = PDGLibrary::Instance()->Find(diquark)->Mass();
160 if( W <= m_frag + m_diquark + fMinESinglet ) {
161 LOG("LeptoHad", pWARN) << "Low invariant mass, W = " << W << " GeV! Returning a null list";
162 LOG("LeptoHad", pWARN) << "frag_quark = " << frag_quark << " -> m = " << m_frag;
163 LOG("LeptoHad", pWARN) << "diquark = " << diquark << " -> m = " << m_diquark;
164 return 0;
165 }
166 // Input the two particles to PYTHIA back to back in the CM frame
167 double e_frag = (W*W - m_diquark*m_diquark + m_frag*m_frag)/2./W;
168 double e_diquark = (W*W + m_diquark*m_diquark - m_frag*m_frag)/2./W;
169 fPythia->Py1ent( -1, frag_quark, e_frag, 0., 0. ); //k(1,2) = 2
170 fPythia->Py1ent( fPythia->GetN()+1, diquark, e_diquark, fPythia->GetPARU(1), 0. ); //k(2,2) = 1
171 }
172
173
174 // If the hit quark is not u or d then is more complicated.
175 // We are using the same procedure use in LEPTO (see lqev.F)
176 // Our initial systemt will look like this -> qqq + hit_q(->frag_q) + rema_q
177 // And we have to input PYTHIA something like this -> frag_q + rema + hadron
178 // These are the posible combinations -> frag_q[q] + meson [qqb] + diquark [qq]
179 // -> frag_q[qb] + baryon [qqq] + quark [q]
180 else {
181
182 // Remnant of the hit quark (which is from the sea) will be of opposite charge
183 int rema_hit_quark = -hit_quark;
184
185 // Check that the trasnfered energy is higher than the mass of the produce quarks plus remnant quark and nucleon
186 double m_frag = PDGLibrary::Instance()->Find(frag_quark)->Mass();
187 double m_rema_hit = PDGLibrary::Instance()->Find(rema_hit_quark)->Mass();
188 if (W <= m_frag + m_rema_hit + 0.9 + fMinESinglet ) {
189 LOG("LeptoHad", pWARN) << "Low invariant mass, W = " << W << " GeV! Returning a null list";
190 LOG("LeptoHad", pWARN) << " frag_quark = " << frag_quark << " -> m = " << m_frag;
191 LOG("LeptoHad", pWARN) << " rema_hit_quark = " << rema_hit_quark << " -> m = " << m_rema_hit;
192 return 0;
193 }
194
195 //PDG of the two hadronic particles for the final state
196 int hadron = 0;
197 int rema = 0;
198
199 int ntwoq = isp ? 2 : 1; //proton two ups & neutron one up
200 int counter = 0;
201
202 // Here we select the id and kinematics of the hadron and rema particles
203 // Some combinations can be kinematically forbiden so we repeat this process
204 // up to 100 times before the event is discarded.
205 while( counter<fMaxIterHad ) {
206
207 // Loop to create a combination of hadron + rema. Two options are possible:
208 // 1) diquark [qq] + meson [qqb]
209 // 2) quark [q] + baryon [qqq]
210 while(hadron==0) {
211 //choose a valence quark and the remaining will be a diquark system
212 int valquark = int(1.+ntwoq/3.+rnd->RndHadro().Rndm());
213 int diquark = 0;
214 if ( valquark==ntwoq ) diquark = rnd->RndHadro().Rndm()>0.75 ? kPdgUDDiquarkS1 : kPdgUDDiquarkS0;
215 else diquark = 1000*ntwoq+100*ntwoq+3;
216
217 // Choose flavours using PYTHIA tool
218 int idum;
219 if ( rema_hit_quark>0 ) { //create a baryon (qqq)
220 pykfdi_(&diquark,&rema_hit_quark,&idum,&hadron);
221 rema = valquark;
222 }
223 else { //create a meson (qqbar)
224 pykfdi_(&valquark,&rema_hit_quark,&idum,&hadron);
225 rema = diquark;
226 }
227 }
228
229 double m_hadron = PDGLibrary::Instance()->Find(hadron)->Mass();
230 double m_rema = PDGLibrary::Instance()->Find(rema)->Mass();
231
232 // Give balancing pT to hadron and rema particles
233 double pT = fRemnantPT * TMath::Sqrt( -1*TMath::Log( rnd->RndHadro().Rndm() ) );
234 double pT2 = TMath::Power(pT,2);
235 double pr = TMath::Power(m_hadron,2)+pT2;
236 //to generate the longitudinal scaling variable z in jet fragmentation using PYTHIA function
237 // Split energy-momentum of remnant using PYTHIA function
238 // z=E-pz fraction for rema forming jet-system with frag_q
239 // 1-z=E-pz fraction for hadron
240 double z;
241 int kfl1 = 1;
242 int kfl3 = 0;
243 pyzdis_(&kfl1,&kfl3,&pr,&z);
244
245 // Energy of trasnfered to the hadron
246 double tm_hadron = pr / z / W;
247 double E_hadron = 0.5 * ( z*W + tm_hadron ); //E_hadron - pz = zW
248 double E_pz = W - tm_hadron;
249 double WT = (1-z) * W * E_pz - pT2;
250
251 // Check if energy in jet system is enough for fragmentation.
252 if ( WT > TMath::Power(m_frag+m_rema+fMinESinglet,2) ) {
253
254 // Energy of transfered to the fragmented quark and rema system
255 // Applying energy conservation
256 WT = TMath::Sqrt( WT + pT2 );
257 double tm_rema = TMath::Power(m_rema,2) + pT2;
258 double E_frag = 0.5 * ( WT + ( TMath::Power(m_frag,2) - tm_rema)/WT ); //E_frag + E_rema = WT
259 double E_rema = 0.5 * ( WT + (-TMath::Power(m_frag,2) + tm_rema)/WT );
260 double x_rema = -1 * TMath::Sqrt( TMath::Power(E_rema,2) - tm_rema );
261 double theta_rema;
262 theta_rema = pyangl_(&x_rema,&pT);
263
264 // Select a phi angle between between particles randomly
265 double phi = 2*kPi*rnd->RndHadro().Rndm();
266
267 double dbez = (E_pz-(1-z)*W)/(E_pz+(1-z)*W);
268 double pz_hadron = -0.5 * ( z*W - tm_hadron );
269
270 // Input the three particles to PYTHIA in the CM frame
271 // If a top quark is produced we decay it because it does not hadronize
272
273 fPythia->Py1ent( -1, frag_quark, E_frag, 0., 0. ); //k(1,2) = 2
274 if (TMath::Abs(frag_quark) > 5 ) {
275 int ip = 1;
276 pydecy_(&ip);
277 }
278 fPythia->Py1ent( fPythia->GetN()+1, rema, E_rema, theta_rema, phi ); //k(2,2) = 1
279
280 int imin = 0;
281 int imax = 0;
282 double the = 0.; double ph = 0.;
283 double dbex = 0.; double dbey = 0.;
284 pyrobo_( &imin , &imax, &the, &ph, &dbex, &dbey , &dbez );
285 double theta_hadron = pyangl_(&pz_hadron,&pT);
286
287 fPythia->SetMSTU( 10, 1 ); //keep the mass value stored in P(I,5), whatever it is.
288 fPythia->SetP( fPythia->GetN()+1, 5, m_hadron );
289 fPythia->Py1ent( fPythia->GetN()+1, hadron, E_hadron, theta_hadron, phi + kPi );
290 fPythia->SetMSTU( 10, 2 ); //find masses according to mass tables as usual.
291
292 // Target remnants required to go backwards in hadronic cms
293 if ( fPythia->GetP(fPythia->GetN()-1,3)<0 && fPythia->GetP(fPythia->GetN(),3)<0 ) break; //quit the while from line 368
294
295 LOG("LeptoHad", pINFO) << "Not backward hadron or rema";
296 LOG("LeptoHad", pINFO) << "hadron = " << hadron << " -> Pz = " << fPythia->GetP(fPythia->GetN(),3) ;
297 LOG("LeptoHad", pINFO) << "rema = " << rema << " -> Pz = " << fPythia->GetP(fPythia->GetN()-1,3) ;
298
299 }
300 else {
301 LOG("LeptoHad", pINFO) << "Low WT value ... ";
302 LOG("LeptoHad", pINFO) << "WT = " << TMath::Sqrt(WT) << " // m_frag = " << m_frag << " // m_rema = " << m_rema;
303 }
304
305 LOG("LeptoHad", pINFO) << "Hadronization paricles not suitable. Trying again... " << counter;
306 counter++;
307 if (counter==100) {
308 LOG("LeptoHad", pWARN) << "Hadronization particles failed after " << counter << " iterations! Returning a null list";
309 return 0;
310 }
311
312 }
313
314 }
315
316 // Introduce a primordial kT system
317 double pT = fPrimordialKT * TMath::Sqrt( -1*TMath::Log( rnd->RndHadro().Rndm() ) );
318 double phi = -2*kPi*rnd->RndHadro().Rndm();
319 double theta = 0.;
320
321 int imin = 0;
322 int imax = 0;
323 double dbex = 0.; double dbey = 0.; double dbez = 0;
324 pyrobo_( &imin , &imax, &theta, &phi, &dbex, &dbey , &dbez );
325 phi = -1 * phi;
326 theta = TMath::ATan(2.*pT/W);
327 pyrobo_( &imin , &imax, &theta, &phi, &dbex, &dbey , &dbez );
328
329 // Run PYTHIA with the input particles
330 fPythia->Pyexec();
331 // Use for debugging purposes
332 //fPythia->Pylist(3);
333 fPythia->GetPrimaries();
334 TClonesArray * pythia_particles = (TClonesArray *) fPythia->ImportParticles("All");
335 // copy PYTHIA container to a new TClonesArray so as to transfer ownership
336 // of the container and of its elements to the calling method
337 int np = pythia_particles->GetEntries();
338 assert(np>0);
339 TClonesArray * particle_list = new TClonesArray("genie::GHepParticle", np);
340 particle_list->SetOwner(true);
341
342 // Boost velocity HCM -> LAB
343 long double beta = p4Hadlong.P()/p4Hadlong.E();
344
345 //fix numbering for events with top
346 bool isTop = false;
347
348 //-- Translate the fragmentation products from TMCParticles to
349 // GHepParticles and copy them to the event record.
350 int mom = event->FinalStateHadronicSystemPosition();
351 assert(mom!=-1);
352
353 TMCParticle * p = 0;
354 TIter particle_iter(pythia_particles);
355 while( (p = (TMCParticle *) particle_iter.Next()) ) {
356
357 int pdgc = p->GetKF();
358 int ks = p->GetKS();
359
360 // Final state particles can not be quarks or diquarks but colorless
361 if(ks == 1) {
362 if( pdg::IsQuark(pdgc) || pdg::IsDiQuark(pdgc) ) {
363 LOG("LeptoHad", pERROR) << "Hadronization failed! Bare quark/di-quarks appear in final state!";
364 return false;
365 }
366 }
367
368 // When top quark is produced, it is immidiately decay before hadronization. Then the decayed
369 // products are hadronized with the hadron remnants. Therefore, we remove the top quark from
370 // the list of particles so that the mother/daugher assigments is at the same level for decayed
371 // products and hadron remnants.
372 if ( pdg::IsTQuark( TMath::Abs(pdgc) ) ) { isTop=true; continue; }
373
374 // fix numbering scheme used for mother/daughter assignments
375 if ( isTop ) {
376 (p->GetParent()==0) ? p->SetParent(p->GetParent() - 1) : p->SetParent(p->GetParent() - 2);
377 p->SetFirstChild (p->GetFirstChild() - 2);
378 p->SetLastChild (p->GetLastChild() - 2);
379 }
380 else {
381 p->SetParent(p->GetParent() - 1);
382 p->SetFirstChild (p->GetFirstChild() - 1);
383 p->SetLastChild (p->GetLastChild() - 1);
384 }
385
386 LongLorentzVector p4long( p->GetPx(), p->GetPy(), p->GetPz(), p->GetEnergy() );
387 p4long.BoostZ(beta);
388 p4long.Rotate(p4Hadlong);
389
390 // Translate from long double to double
391 TLorentzVector p4( (double)p4long.Px(), (double)p4long.Py(), (double)p4long.Pz(), (double)p4long.E() );
392
393 // Somtimes PYTHIA output particles with E smaller than its mass. This is wrong,
394 // so we assume that the are at rest.
395 double massPDG = PDGLibrary::Instance()->Find(pdgc)->Mass();
396 if ( (ks==1 || ks==4) && p4.E()<massPDG ) {
397 LOG("LeptoHad", pINFO) << "Putting at rest one stable particle generated by PYTHIA because E < m";
398 LOG("LeptoHad", pINFO) << "PDG = " << pdgc << " // State = " << ks;
399 LOG("LeptoHad", pINFO) << "E = " << p4.E() << " // |p| = " << p4.P();
400 LOG("LeptoHad", pINFO) << "p = [ " << p4.Px() << " , " << p4.Py() << " , " << p4.Pz() << " ]";
401 LOG("LeptoHad", pINFO) << "m = " << p4.M() << " // mpdg = " << massPDG;
402 p4.SetXYZT(0,0,0,massPDG);
403 }
404
405 // copy final state particles to the event record
407
408 int im = mom + 1 + p->GetParent();
409 int ifc = (p->GetFirstChild() <= -1) ? -1 : mom + 1 + p->GetFirstChild();
410 int ilc = (p->GetLastChild() <= -1) ? -1 : mom + 1 + p->GetLastChild();
411
412 double vx = vtx.X() + p->GetVx()*1e12; //pythia gives position in [mm] while genie uses [fm]
413 double vy = vtx.Y() + p->GetVy()*1e12;
414 double vz = vtx.Z() + p->GetVz()*1e12;
415 double vt = vtx.T() + p->GetTime()*(units::millimeter/units::second);
416 TLorentzVector pos( vx, vy, vz, vt );
417
418 event->AddParticle( pdgc, ist, im,-1, ifc, ilc, p4, pos );
419
420 }
421
422 return true;
423#else
424 return false;
425#endif
426
427}
#define pINFO
Definition Messenger.h:62
#define pERROR
Definition Messenger.h:59
#define pDEBUG
Definition Messenger.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
#define pWARN
Definition Messenger.h:60
const TLorentzVector * P4(void) const
virtual GHepParticle * Probe(void) const
virtual GHepParticle * FinalStatePrimaryLepton(void) const
virtual int HitNucleonPosition(void) const
virtual GHepParticle * HitNucleon(void) const
const Target & Tgt(void) const
const XclsTag & ExclTag(void) const
Definition Interaction.h:72
const Kinematics & Kine(void) const
Definition Interaction.h:71
const InitialState & InitState(void) const
Definition Interaction.h:69
Kinematics * KinePtr(void) const
Definition Interaction.h:76
void SetHadSystP4(const TLorentzVector &p4)
double W(bool selected=false) const
void SetW(double W, bool selected=false)
static PDGLibrary * Instance(void)
TParticlePDG * Find(int pdgc, bool must_exist=true)
TRandom3 & RndHadro(void) const
rnd number generator used by hadronization models
Definition RandomGen.h:53
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
int HitNucPdg(void) const
Definition Target.cxx:304
int HitQrkPdg(void) const
Definition Target.cxx:242
bool HitQrkIsSet(void) const
Definition Target.cxx:292
int FinalQuarkPdg(void) const
Definition XclsTag.h:72
bool IsQuark(int pdgc)
Definition PDGUtils.cxx:250
bool IsProton(int pdgc)
Definition PDGUtils.cxx:336
bool IsTQuark(int pdgc)
Definition PDGUtils.cxx:291
bool IsUQuark(int pdgc)
Definition PDGUtils.cxx:266
bool IsDQuark(int pdgc)
Definition PDGUtils.cxx:271
bool IsDiQuark(int pdgc)
Definition PDGUtils.cxx:231
static constexpr double millimeter
Definition Units.h:41
static constexpr double second
Definition Units.h:37
double W(const Interaction *const i)
@ kIStDISPreFragmHadronicState
Definition GHepStatus.h:35
@ kIStStableFinalState
Definition GHepStatus.h:30
const int kPdgUDDiquarkS0
Definition PDGCodes.h:56
const int kPdgUUDiquarkS1
Definition PDGCodes.h:58
const int kPdgHadronicSyst
Definition PDGCodes.h:210
enum genie::EGHepStatus GHepStatus_t
const int kPdgUDDiquarkS1
Definition PDGCodes.h:57
const int kPdgDDDiquarkS1
Definition PDGCodes.h:55

References genie::utils::math::LongLorentzVector::BoostZ(), genie::utils::math::LongLorentzVector::Dx(), genie::utils::math::LongLorentzVector::Dy(), genie::utils::math::LongLorentzVector::Dz(), genie::utils::math::LongLorentzVector::E(), genie::Interaction::ExclTag(), genie::XclsTag::FinalQuarkPdg(), genie::PDGLibrary::Find(), fMaxIterHad, fMinESinglet, fPrimordialKT, fRemnantPT, fWmin, genie::Target::HitNucPdg(), genie::Target::HitQrkIsSet(), genie::Target::HitQrkPdg(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::RandomGen::Instance(), genie::pdg::IsDiQuark(), genie::pdg::IsDQuark(), genie::pdg::IsProton(), genie::pdg::IsQuark(), genie::pdg::IsTQuark(), genie::pdg::IsUQuark(), genie::Interaction::Kine(), genie::Interaction::KinePtr(), genie::kIStDISPreFragmHadronicState, genie::kIStStableFinalState, genie::kPdgDDDiquarkS1, genie::kPdgHadronicSyst, genie::kPdgUDDiquarkS0, genie::kPdgUDDiquarkS1, genie::kPdgUUDiquarkS1, genie::constants::kPi, LOG, genie::utils::math::LongLorentzVector::M(), genie::utils::math::LongLorentzVector::M2(), genie::units::millimeter, genie::utils::math::LongLorentzVector::P(), pDEBUG, pERROR, pINFO, pWARN, genie::utils::math::LongLorentzVector::Px(), genie::utils::math::LongLorentzVector::Py(), genie::utils::math::LongLorentzVector::Pz(), genie::RandomGen::RndHadro(), genie::utils::math::LongLorentzVector::Rotate(), genie::units::second, genie::Kinematics::SetHadSystP4(), genie::Kinematics::SetW(), genie::InitialState::Tgt(), and genie::Kinematics::W().

Referenced by ProcessEventRecord().

◆ Initialize()

void LeptoHadPythia6::Initialize ( void ) const
private

Definition at line 499 of file LeptoHadPythia6.cxx.

500{
501#ifdef __GENIE_PYTHIA6_ENABLED__
502 fPythia = TPythia6::Instance();
503
504 // sync GENIE/PYTHIA6 seed number
506#endif
507
508}

References genie::RandomGen::Instance().

Referenced by LeptoHadPythia6(), and LeptoHadPythia6().

◆ LoadConfig()

void LeptoHadPythia6::LoadConfig ( void )
private

Definition at line 441 of file LeptoHadPythia6.cxx.

442{
443
444#ifdef __GENIE_PYTHIA6_ENABLED__
445 GetParam("MaxIter-Had", fMaxIterHad ) ;
446
447 // Width of Gaussian distribution for transverse momentums
448 // Define in LEPTO with PARL(3) and PARL(14)
449 GetParam("Primordial-kT", fPrimordialKT ) ;
450 GetParam("Remnant-pT", fRemnantPT ) ;
451
452 // It is, with quark masses added, used to define the minimum allowable energy of a colour-singlet parton system.
453 GetParam( "Energy-Singlet", fMinESinglet ) ;
454
455 GetParam( "Xsec-Wmin", fWmin ) ;
456
457 GetParam( "SSBarSuppression", fSSBarSuppression );
458 GetParam( "GaussianPt2", fGaussianPt2 );
459 GetParam( "NonGaussianPt2Tail", fNonGaussianPt2Tail );
460 GetParam( "RemainingEnergyCutoff", fRemainingECutoff );
461 GetParam( "DiQuarkSuppression", fDiQuarkSuppression );
462 GetParam( "LightVMesonSuppression", fLightVMesonSuppression );
463 GetParam( "SVMesonSuppression", fSVMesonSuppression );
464 GetParam( "Lunda", fLunda );
465 GetParam( "Lundb", fLundb );
466 GetParam( "LundaDiq", fLundaDiq );
467
468 int warnings; GetParam( "Warnings", warnings ) ;
469 int errors; GetParam( "Errors", errors ) ;
470 int qrk_mass; GetParam( "QuarkMass", qrk_mass ) ;
471
472 // PYTHIA6 parameters only valid for HEDIS
473 fPythia->SetPARP(2, fWmin); // (D = 10. GeV) lowest c.m. energy for the event as a whole that the program will accept to simulate. (bellow 2GeV pythia crashes)
474 fPythia->SetMSTU(26, warnings); // (Default=10) maximum number of warnings that are printed
475 fPythia->SetMSTU(22, errors); // (Default=10) maximum number of errors that are printed
476 fPythia->SetMSTJ(93, qrk_mass); // light (d, u, s, c, b) quark masses are taken from PARF(101) - PARF(105) rather than PMAS(1,1) - PMAS(5,1). Diquark masses are given as sum of quark masses, without spin splitting term.
477 fPythia->SetPMAS(24,1,kMw); // mass of the W boson (pythia=80.450 // genie=80.385)
478 fPythia->SetPMAS(24,2,0.); // set to 0 the width of the W boson to avoid problems with energy conservation
479 fPythia->SetPMAS(6,2,0.); // set to 0 the width of the top to avoid problems with energy conservation
480 fPythia->SetMDME(192,1,0); // W->dbar+t decay off
481 fPythia->SetMDME(196,1,0); // W->cbar+t decay off
482 fPythia->SetMDME(200,1,0); // W->cbar+t decay off
483
484 // PYTHIA tuned parameters
485 fPythia->SetPARJ(2, fSSBarSuppression );
486 fPythia->SetPARJ(21, fGaussianPt2 );
487 fPythia->SetPARJ(23, fNonGaussianPt2Tail );
488 fPythia->SetPARJ(33, fRemainingECutoff );
489 fPythia->SetPARJ(1, fDiQuarkSuppression );
490 fPythia->SetPARJ(11, fLightVMesonSuppression );
491 fPythia->SetPARJ(12, fSVMesonSuppression );
492 fPythia->SetPARJ(41, fLunda );
493 fPythia->SetPARJ(42, fLundb );
494 fPythia->SetPARJ(45, fLundaDiq );
495#endif
496
497}
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
double fDiQuarkSuppression
di-quark suppression parameter
double fGaussianPt2
gaussian pt2 distribution width
double fLunda
Lund a parameter.
double fLundaDiq
adjustment of Lund a for di-quark
double fNonGaussianPt2Tail
non gaussian pt2 tail parameterization
double fSVMesonSuppression
strange vector meson suppression
double fSSBarSuppression
ssbar suppression
double fLundb
Lund b parameter.
double fRemainingECutoff
remaining E cutoff stopping fragmentation
double fLightVMesonSuppression
light vector meson suppression

References fDiQuarkSuppression, fGaussianPt2, fLightVMesonSuppression, fLunda, fLundaDiq, fLundb, fMaxIterHad, fMinESinglet, fNonGaussianPt2Tail, fPrimordialKT, fRemainingECutoff, fRemnantPT, fSSBarSuppression, fSVMesonSuppression, fWmin, genie::Algorithm::GetParam(), and genie::constants::kMw.

Referenced by Configure(), and Configure().

◆ ProcessEventRecord()

void LeptoHadPythia6::ProcessEventRecord ( GHepRecord * event) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 48 of file LeptoHadPythia6.cxx.

49{
50
51 if(!this->Hadronize(event)) {
52 LOG("LeptoHad", pWARN) << "Hadronization failed!";
53 event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
54 genie::exceptions::EVGThreadException exception;
55 exception.SetReason("Could not simulate the hadronic system");
56 exception.SwitchOnFastForward();
57 throw exception;
58 return;
59 }
60
61
62}
bool Hadronize(GHepRecord *event) const
@ kHadroSysGenErr
Definition GHepFlags.h:32

References Hadronize(), genie::kHadroSysGenErr, LOG, pWARN, genie::exceptions::EVGThreadException::SetReason(), and genie::exceptions::EVGThreadException::SwitchOnFastForward().

Member Data Documentation

◆ fDiQuarkSuppression

double genie::LeptoHadPythia6::fDiQuarkSuppression
private

di-quark suppression parameter

Definition at line 79 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fGaussianPt2

double genie::LeptoHadPythia6::fGaussianPt2
private

gaussian pt2 distribution width

Definition at line 76 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fLightVMesonSuppression

double genie::LeptoHadPythia6::fLightVMesonSuppression
private

light vector meson suppression

Definition at line 80 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fLunda

double genie::LeptoHadPythia6::fLunda
private

Lund a parameter.

Definition at line 82 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fLundaDiq

double genie::LeptoHadPythia6::fLundaDiq
private

adjustment of Lund a for di-quark

Definition at line 84 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fLundb

double genie::LeptoHadPythia6::fLundb
private

Lund b parameter.

Definition at line 83 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fMaxIterHad

int genie::LeptoHadPythia6::fMaxIterHad
private

Definition at line 68 of file LeptoHadPythia6.h.

Referenced by Hadronize(), and LoadConfig().

◆ fMinESinglet

double genie::LeptoHadPythia6::fMinESinglet
private

Definition at line 71 of file LeptoHadPythia6.h.

Referenced by Hadronize(), and LoadConfig().

◆ fNonGaussianPt2Tail

double genie::LeptoHadPythia6::fNonGaussianPt2Tail
private

non gaussian pt2 tail parameterization

Definition at line 77 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fPrimordialKT

double genie::LeptoHadPythia6::fPrimordialKT
private

Definition at line 69 of file LeptoHadPythia6.h.

Referenced by Hadronize(), and LoadConfig().

◆ fRemainingECutoff

double genie::LeptoHadPythia6::fRemainingECutoff
private

remaining E cutoff stopping fragmentation

Definition at line 78 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fRemnantPT

double genie::LeptoHadPythia6::fRemnantPT
private

Definition at line 70 of file LeptoHadPythia6.h.

Referenced by Hadronize(), and LoadConfig().

◆ fSSBarSuppression

double genie::LeptoHadPythia6::fSSBarSuppression
private

ssbar suppression

Definition at line 75 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fSVMesonSuppression

double genie::LeptoHadPythia6::fSVMesonSuppression
private

strange vector meson suppression

Definition at line 81 of file LeptoHadPythia6.h.

Referenced by LoadConfig().

◆ fWmin

double genie::LeptoHadPythia6::fWmin
private

Definition at line 72 of file LeptoHadPythia6.h.

Referenced by Hadronize(), and LoadConfig().


The documentation for this class was generated from the following files: