GENIEGenerator
Loading...
Searching...
No Matches
Pythia8Hadro2019.cxx
Go to the documentation of this file.
1//____________________________________________________________________________
2/*
3 Copyright (c) 2003-2025, The GENIE Collaboration
4 For the full text of the license visit http://copyright.genie-mc.org
5
6 Costas Andreopoulos <c.andreopoulos \at cern.ch>
7 University of Liverpool
8
9 Shivesh Mandalia <s.p.mandalia@qmul.ac.uk>
10 Queen Mary University of London
11*/
12//____________________________________________________________________________
13
14#include <RVersion.h>
15#include <TClonesArray.h>
16// Avoid the inclusion of dlfcn.h by Pythia.h that CINT is not able to process
17#ifdef __CINT__
18#define _DLFCN_H_
19#define _DLFCN_H
20#endif
21
35
36#ifdef __GENIE_PYTHIA8_ENABLED__
37#include "Pythia8/Pythia.h"
38#endif
39
40using namespace genie;
41using namespace genie::constants;
42
43//____________________________________________________________________________
45PythiaBaseHadro2019("genie::Pythia8Hadro2019")
46{
47 this->Initialize();
48}
49//____________________________________________________________________________
51PythiaBaseHadro2019("genie::Pythia8Hadro2019", config)
52{
53 this->Initialize();
54}
55//____________________________________________________________________________
60//____________________________________________________________________________
62#ifdef __GENIE_PYTHIA8_ENABLED__
63 event // avoid unused variable warning if PYTHIA6 is not enabled
64#endif
65) const
66{
67#ifdef __GENIE_PYTHIA8_ENABLED__
69#else
70 LOG("Pythia8Had", pFATAL)
71 << "Calling GENIE/PYTHIA8 hadronization modules without enabling PYTHIA8";
72 gAbortingInErr = true;
73 std::exit(1);
74#endif
75}
76//____________________________________________________________________________
78#ifdef __GENIE_PYTHIA8_ENABLED__
79 event // avoid unused variable warning if PYTHIA6 is not enabled
80#endif
81) const
82{
83#ifdef __GENIE_PYTHIA8_ENABLED__
84 LOG("Pythia8Had", pNOTICE) << "Running PYTHIA8 hadronizer";
85
86 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
87
88 const Interaction * interaction = event->Summary();
89 const Kinematics & kinematics = interaction->Kine();
90 double W = kinematics.W();
91
92 LOG("Pythia8Had", pNOTICE)
93 << "Fragmentation: "
94 << "q = " << fLeadingQuark << ", qq = " << fRemnantDiquark
95 << ", W = " << W << " GeV";
96
97 // Hadronize
98
99 LOG("Pythia8Had", pDEBUG) << "Reseting PYTHIA8 event";
100 gPythia->event.reset();
101
102 // Get quark/diquark masses
103 double mA = gPythia->particleData.m0(fLeadingQuark);
104 double mB = gPythia->particleData.m0(fRemnantDiquark);
105
106 LOG("Pythia8Had", pINFO)
107 << "Leading quark mass = " << mA
108 << " GeV, remnant diqurak mass = " << mB << ", GeV";
109
110 // Calculate quark/diquark energy/momentum
111 double pzAcm = 0.5 * Pythia8::sqrtpos(
112 (W + mA + mB) * (W - mA - mB) * (W - mA + mB) * (W + mA - mB) ) / W;
113 double pzBcm = -pzAcm;
114 double eA = sqrt(mA*mA + pzAcm*pzAcm);
115 double eB = sqrt(mB*mB + pzBcm*pzBcm);
116
117 LOG("Pythia8Had", pINFO)
118 << "Quark: (pz = " << pzAcm << ", E = " << eA << ") GeV, "
119 << "Diquark: (pz = " << pzBcm << ", E = " << eB << ") GeV";
120
121 // Pythia8 status code for outgoing particles of the hardest subprocesses is 23
122 // anti/colour tags for these 2 particles must complement each other
123 LOG("Pythia8Had", pDEBUG) << "Appending quark/diquark into the PYTHIA8 event";
124 gPythia->event.append(fLeadingQuark, 23, 101, 0, 0., 0., pzAcm, eA, mA);
125 gPythia->event.append(fRemnantDiquark, 23, 0, 101, 0., 0., pzBcm, eB, mB);
126 gPythia->event.list();
127
128 LOG("Pythia8Had", pDEBUG) << "Generating next PYTHIA8 event";
129 gPythia->next();
130
131 // List the event information
132 //gPythia->event.list();
133 //gPythia->stat();
134
135 // Get LUJETS record
136 LOG("Pythia8Had", pDEBUG) << "Copying PYTHIA8 event record into GENIE's";
137 Pythia8::Event &fEvent = gPythia->event;
138 int np = fEvent.size();
139 assert(np>0);
140
141 TClonesArray * particle_list = new TClonesArray("genie::GHepParticle", np);
142 particle_list->SetOwner(true);
143
144 // Hadronic 4vec
145 TLorentzVector p4Had = kinematics.HadSystP4();
146
147 // Vector defining rotation from LAB to LAB' (z:= \vec{phad})
148 TVector3 unitvq = p4Had.Vect().Unit();
149
150 // Boost velocity LAB' -> HCM
151 TVector3 beta(0,0,p4Had.P()/p4Had.Energy());
152
153 // Check target and decide appropriate status code for f/s particles
154 bool is_nucleus = interaction->InitState().Tgt().IsNucleus();
156
157 // Get the index of the mother of the hadronic system
158 int mom = event->FinalStateHadronicSystemPosition();
159 assert(mom!=-1);
160
161 // Get the neutrino vertex position (all hadrons positions set to this point)
162 GHepParticle * neutrino = event->Probe();
163 const TLorentzVector & vtx = *(neutrino->X4());
164
165 // Loop over PYTHIA8 event particles and copy relevant entries
166 for (int i = 0; i < np; i++) {
167
168 if (fEvent[i].id() == 90) continue; // ignore (system) pseudoparticle
169
170 // Get PYTHIA8 particle ID and status codes
171 int particle_pdg_code = fEvent[i].id();
172 int pythia_particle_status = fEvent[i].status();
173
174 // Sanity check
175 if(pythia_particle_status > 0) {
176 if( pdg::IsQuark (particle_pdg_code) ||
177 pdg::IsDiQuark(particle_pdg_code) )
178 {
179 LOG("Pythia8Had", pERROR)
180 << "Hadronization failed! Bare quarks appear in final state!";
181 return false;
182 }
183 }
184
185 // Copy the initial q, qq (status = -23) and all undecayed particles
186 // Ignore particles we asked PYTHIA to decay (eg omega, Delta) but
187 // record their decay products
188 bool copy = (pythia_particle_status==-23) || (pythia_particle_status > 0);
189 if(copy) {
190 // The fragmentation products are generated in the hadronic CM frame
191 // where the z>0 axis is the \vec{phad} direction. For each particle
192 // returned by the hadronizer:
193 // - boost it back to LAB' frame {z:=\vec{phad}} / doesn't affect pT
194 // - rotate its 3-momentum from LAB' to LAB
195 TLorentzVector p4o(
196 fEvent[i].px(), fEvent[i].py(), fEvent[i].pz(), fEvent[i].e());
197 p4o.Boost(beta);
198 TVector3 p3 = p4o.Vect();
199 p3.RotateUz(unitvq);
200 TLorentzVector p4(p3,p4o.Energy());
201
202 // Set the proper GENIE status according to a number of things:
203 // interaction on a nucleus or nucleon, particle type
204 GHepStatus_t ist = (pythia_particle_status > 0) ?
206 // Handle gammas, and leptons that might come from internal pythia decays
207 // mark them as final state particles
208 bool is_gamma = (particle_pdg_code == kPdgGamma);
209 bool is_nu = pdg::IsNeutralLepton(particle_pdg_code);
210 bool is_lchg = pdg::IsChargedLepton(particle_pdg_code);
211 bool not_hadr = is_gamma || is_nu || is_lchg;
212 if(not_hadr) { ist = kIStStableFinalState; }
213
214 // Set mother/daugher indices
215 // int mother1 = mom+ fEvent[i].mother1();
216 // int mother2 = (pythia_particle_status > 0) ? mom + fEvent[i].mother2() : -1;
217 int mother1 = mom; // fEvent[i].mother1();
218 int mother2 = -1; //(pythia_particle_status > 0) ? mom + fEvent[i].mother2() : -1;
219 if(pythia_particle_status > 0) {
220 mother1 = mom+1;
221 mother2 = mom+2;
222 }
223 int daughter1 = -1;//(fEvent[i].daughter1() <= 0 ) ? -1 : mom + fEvent[i].daughter1();
224 int daughter2 = -1;//(fEvent[i].daughter1() <= 0 ) ? -1 : mom + fEvent[i].daughter2();
225
226 // Create GHepParticle
227 GHepParticle particle = GHepParticle(
228 particle_pdg_code, // pdg
229 ist, // status
230 mother1, // first parent
231 mother2, // second parent
232 daughter1, // first daughter
233 daughter2, // second daughter
234 p4.Px(), // px
235 p4.Py(), // py
236 p4.Pz(), // pz
237 p4.Energy(), // e
238 vtx.X(), // x
239 vtx.Y(), // y
240 vtx.Z(), // z
241 vtx.T() // t
242 );
243
244 LOG("Pythia8Had", pDEBUG)
245 << "Adding final state particle pdgc = " << particle.Pdg()
246 << " with status = " << particle.Status();
247
248 // Insert the particle in the list
249 event->AddParticle(particle);
250 }// copy?
251 }// loop over particles
252
253 return true;
254
255#else
256 return false;
257#endif
258}
259//____________________________________________________________________________
261{
262#ifdef __GENIE_PYTHIA8_ENABLED__
263 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
264
265 fOriDecayFlag_pi0 = gPythia->particleData.canDecay(kPdgPi0);
266 fOriDecayFlag_K0 = gPythia->particleData.canDecay(kPdgK0);
267 fOriDecayFlag_K0b = gPythia->particleData.canDecay(kPdgAntiK0);
268 fOriDecayFlag_L0 = gPythia->particleData.canDecay(kPdgLambda);
269 fOriDecayFlag_L0b = gPythia->particleData.canDecay(kPdgAntiLambda);
270 fOriDecayFlag_Dm = gPythia->particleData.canDecay(kPdgP33m1232_DeltaM);
271 fOriDecayFlag_D0 = gPythia->particleData.canDecay(kPdgP33m1232_Delta0);
272 fOriDecayFlag_Dp = gPythia->particleData.canDecay(kPdgP33m1232_DeltaP);
273 fOriDecayFlag_Dpp = gPythia->particleData.canDecay(kPdgP33m1232_DeltaPP);
274
275#ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
276 LOG("Pythia8Had", pDEBUG)
277 << "Original PYTHIA6 decay flags:"
278 << "\n pi0 = " << fOriDecayFlag_pi0
279 << "\n K0 = " << fOriDecayFlag_K0
280 << "\n \bar{K0} = " << fOriDecayFlag_K0b
281 << "\n Lambda = " << fOriDecayFlag_L0
282 << "\n \bar{Lambda0} = " << fOriDecayFlag_L0b
283 << "\n D- = " << fOriDecayFlag_Dm
284 << "\n D0 = " << fOriDecayFlag_D0
285 << "\n D+ = " << fOriDecayFlag_Dp
286 << "\n D++ = " << fOriDecayFlag_Dpp;
287#endif
288
289#endif
290}
291//____________________________________________________________________________
293{
294#ifdef __GENIE_PYTHIA8_ENABLED__
295 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
296
297 gPythia->particleData.mayDecay(kPdgPi0, fReqDecayFlag_pi0 );
298 gPythia->particleData.mayDecay(kPdgK0, fReqDecayFlag_K0 );
299 gPythia->particleData.mayDecay(kPdgAntiK0, fReqDecayFlag_K0b );
300 gPythia->particleData.mayDecay(kPdgLambda, fReqDecayFlag_L0 );
301 gPythia->particleData.mayDecay(kPdgAntiLambda, fReqDecayFlag_L0b );
302 gPythia->particleData.mayDecay(kPdgP33m1232_DeltaM, fReqDecayFlag_Dm );
303 gPythia->particleData.mayDecay(kPdgP33m1232_Delta0, fReqDecayFlag_D0 );
304 gPythia->particleData.mayDecay(kPdgP33m1232_DeltaP, fReqDecayFlag_Dp );
305 gPythia->particleData.mayDecay(kPdgP33m1232_DeltaPP, fReqDecayFlag_Dpp );
306#endif
307}
308//____________________________________________________________________________
310{
311#ifdef __GENIE_PYTHIA8_ENABLED__
312 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
313
314 gPythia->particleData.mayDecay(kPdgPi0, fOriDecayFlag_pi0 );
315 gPythia->particleData.mayDecay(kPdgK0, fOriDecayFlag_K0 );
316 gPythia->particleData.mayDecay(kPdgAntiK0, fOriDecayFlag_K0b );
317 gPythia->particleData.mayDecay(kPdgLambda, fOriDecayFlag_L0 );
318 gPythia->particleData.mayDecay(kPdgAntiLambda, fOriDecayFlag_L0b );
319 gPythia->particleData.mayDecay(kPdgP33m1232_DeltaM, fOriDecayFlag_Dm );
320 gPythia->particleData.mayDecay(kPdgP33m1232_Delta0, fOriDecayFlag_D0 );
321 gPythia->particleData.mayDecay(kPdgP33m1232_DeltaP, fOriDecayFlag_Dp );
322 gPythia->particleData.mayDecay(kPdgP33m1232_DeltaPP, fOriDecayFlag_Dpp );
323#endif
324}
325//____________________________________________________________________________
327{
328 Algorithm::Configure(config);
329 this->LoadConfig();
330}
331//____________________________________________________________________________
333{
334 Algorithm::Configure(config);
335 this->LoadConfig();
336}
337//____________________________________________________________________________
339{
341
342#ifdef __GENIE_PYTHIA8_ENABLED__
343 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
344
345 gPythia->settings.parm("StringFlav:probStoUD", fSSBarSuppression);
346 gPythia->settings.parm("Diffraction:primKTwidth", fGaussianPt2);
347 gPythia->settings.parm("StringPT:enhancedFraction", fNonGaussianPt2Tail);
348 gPythia->settings.parm("StringFragmentation:stopMass", fRemainingECutoff);
349 gPythia->settings.parm("StringFlav:probQQtoQ", fDiQuarkSuppression);
350 gPythia->settings.parm("StringFlav:mesonUDvector", fLightVMesonSuppression);
351 gPythia->settings.parm("StringFlav:mesonSvector", fSVMesonSuppression);
352 gPythia->settings.parm("StringZ:aLund", fLunda);
353 gPythia->settings.parm("StringZ:bLund", fLundb);
354 gPythia->settings.parm("StringZ:aExtraDiquark", fLundaDiq);
355
356 gPythia->init(); // needed again to read the above?
357
358#endif
359
360 LOG("Pythia8Had", pDEBUG) << this->GetConfig();
361}
362//____________________________________________________________________________
364{
365#ifdef __GENIE_PYTHIA8_ENABLED__
366 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
367
368 gPythia->readString("ProcessLevel:all = off");
369 gPythia->readString("Print:quiet = on");
370
371 // sync GENIE and PYTHIA8 seeds
373 long int seed = rnd->GetSeed();
374 gPythia->readString("Random:setSeed = on");
375 gPythia->settings.mode("Random:seed", seed);
376 LOG("Pythia8Had", pINFO)
377 << "PYTHIA8 seed = " << gPythia->settings.mode("Random:seed");
378
379 gPythia->init();
380
381#endif
382}
383//____________________________________________________________________________
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#define pERROR
Definition Messenger.h:59
#define pFATAL
Definition Messenger.h:56
#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
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
virtual const Registry & GetConfig(void) const
virtual void Configure(const Registry &config)
Definition Algorithm.cxx:62
STDHEP-like event record entry that can fit a particle or a nucleus.
int Pdg(void) const
const TLorentzVector * X4(void) const
GHepStatus_t Status(void) const
GENIE's GHEP MC event record.
Definition GHepRecord.h:45
const Target & Tgt(void) const
Summary information for an interaction.
Definition Interaction.h:56
const Kinematics & Kine(void) const
Definition Interaction.h:71
const InitialState & InitState(void) const
Definition Interaction.h:69
Generated/set kinematical variables for an event.
Definition Kinematics.h:39
const TLorentzVector & HadSystP4(void) const
Definition Kinematics.h:66
double W(bool selected=false) const
void CopyOriginalDecayFlags(void) const
void ProcessEventRecord(GHepRecord *event) const
bool Hadronize(GHepRecord *event) const
void SetDesiredDecayFlags(void) const
void RestoreOriginalDecayFlags(void) const
void Configure(const Registry &config)
static Pythia8Singleton * Instance(void)
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
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
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition RandomGen.h:29
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
long int GetSeed(void) const
Definition RandomGen.h:82
A registry. Provides the container for algorithm configuration parameters.
Definition Registry.h:65
bool IsNucleus(void) const
Definition Target.cxx:272
const double e
Basic constants.
bool IsQuark(int pdgc)
Definition PDGUtils.cxx:250
bool IsChargedLepton(int pdgc)
Definition PDGUtils.cxx:101
bool IsNeutralLepton(int pdgc)
Definition PDGUtils.cxx:95
bool IsDiQuark(int pdgc)
Definition PDGUtils.cxx:231
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
@ kIStHadronInTheNucleus
Definition GHepStatus.h:37
@ kIStDISPreFragmHadronicState
Definition GHepStatus.h:35
@ kIStStableFinalState
Definition GHepStatus.h:30
bool gAbortingInErr
Definition Messenger.cxx:34
const int kPdgP33m1232_Delta0
Definition PDGCodes.h:105
const int kPdgPi0
Definition PDGCodes.h:160
const int kPdgP33m1232_DeltaPP
Definition PDGCodes.h:107
enum genie::EGHepStatus GHepStatus_t
const int kPdgAntiLambda
Definition PDGCodes.h:86
const int kPdgP33m1232_DeltaM
Definition PDGCodes.h:104
const int kPdgLambda
Definition PDGCodes.h:85
const int kPdgP33m1232_DeltaP
Definition PDGCodes.h:106
const int kPdgGamma
Definition PDGCodes.h:189
const int kPdgAntiK0
Definition PDGCodes.h:175
const int kPdgK0
Definition PDGCodes.h:174