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

Interface to PYTHIA particle decayer.
The Pythia8Decayer2023 is a concrete implementation of the Decayer interface. More...

#include <Pythia8Decayer2023.h>

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

Public Member Functions

 Pythia8Decayer2023 ()
 Pythia8Decayer2023 (string config)
virtual ~Pythia8Decayer2023 ()
void ProcessEventRecord (GHepRecord *event) const

Private Member Functions

void Initialize (void) const
bool IsHandled (int pdgc) const
void InhibitDecay (int pdgc, TDecayChannel *ch=0) const
void UnInhibitDecay (int pdgc, TDecayChannel *ch=0) const
bool Decay (int ip, GHepRecord *event) const
double SumOfBranchingRatios (int kc) const
int FindPythiaDecayChannel (int kc, TDecayChannel *ch) const
bool MatchDecayChannels (int ic, TDecayChannel *ch) const

Private Attributes

double fWeight

Additional Inherited Members

Protected Member Functions inherited from genie::Decayer
 Decayer ()
 Decayer (string name)
 Decayer (string name, string config)
virtual void LoadConfig (void)
virtual bool ToBeDecayed (int pdgc, GHepStatus_t ist) const
virtual bool IsUnstable (int pdgc) const
virtual ~Decayer ()
void Configure (const Registry &config)
void Configure (string config)
Protected Member Functions inherited from genie::EventRecordVisitorI
 EventRecordVisitorI ()
 EventRecordVisitorI (string name)
 EventRecordVisitorI (string name, string config)
virtual ~EventRecordVisitorI ()
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.
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.
Static Protected 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 Attributes inherited from genie::Decayer
bool fGenerateWeighted
 generate weighted or unweighted decays?
bool fRunBefHadroTransp
 is invoked before or after FSI?
PDGCodeList fParticlesToDecay
 list of particles to be decayed
PDGCodeList fParticlesNotToDecay
 list of particles for which decay is inhibited
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

Interface to PYTHIA particle decayer.
The Pythia8Decayer2023 is a concrete implementation of the Decayer interface.

Author
Robert Hatcher <rhatcher \at fnal.gov> Fermilab
Created:\n December 21, 2023
License:\n Copyright (c) 2003-2025, The GENIE Collaboration
For the full text of the license visit http://copyright.genie-mc.org

Definition at line 32 of file Pythia8Decayer2023.h.

Constructor & Destructor Documentation

◆ Pythia8Decayer2023() [1/2]

Pythia8Decayer2023::Pythia8Decayer2023 ( )

Definition at line 37 of file Pythia8Decayer2023.cxx.

37 :
38Decayer("genie::Pythia8Decayer2023")
39{
40 this->Initialize();
41}

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

◆ Pythia8Decayer2023() [2/2]

Pythia8Decayer2023::Pythia8Decayer2023 ( string config)

Definition at line 43 of file Pythia8Decayer2023.cxx.

43 :
44Decayer("genie::Pythia8Decayer2023", config)
45{
46 this->Initialize();
47}

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

◆ ~Pythia8Decayer2023()

Pythia8Decayer2023::~Pythia8Decayer2023 ( )
virtual

Definition at line 49 of file Pythia8Decayer2023.cxx.

50{
51
52}

Member Function Documentation

◆ Decay()

bool Pythia8Decayer2023::Decay ( int ip,
GHepRecord * event ) const
private

Definition at line 87 of file Pythia8Decayer2023.cxx.

88{
89 fWeight = 1.; // reset previous decay weight
90
91#ifdef __GENIE_PYTHIA8_ENABLED__
92
93 // Get the one-and-only instance of Pythia8 that GENIE will use
94 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
95
96 // Get particle to be decayed
97 GHepParticle * decay_particle = event->Particle(decay_particle_id);
98 if(!decay_particle) return 0;
99
100 // Check whether the interaction is off a nuclear target or free nucleon
101 // Depending on whether this module is run before or after the hadron
102 // transport module it would affect the daughters status code
103 GHepParticle * target_nucleus = event->TargetNucleus();
104 bool in_nucleus = (target_nucleus!=0);
105
106 // Get the particle 4-momentum, 4-position and PDG code
107 TLorentzVector decay_particle_p4 = *(decay_particle->P4());
108 TLorentzVector decay_particle_x4 = *(decay_particle->X4());
109 int decay_particle_pdg_code = decay_particle->Pdg();
110
111 TVector3 polz;
112 decay_particle->GetPolarization(polz);
113
114 gPythia->event.reset();
115
116 // check if pdgid is consistent with Pythia8 particle table
117 //
118 if ( !gPythia->particleData.findParticle(decay_particle_pdg_code) ) {
119 LOG("Pythia8Decay", pWARN)
120 << " can NOT find pdgid = " << decay_particle_pdg_code
121 << " in Pythia8::ParticleData";
122 return false;
123 }
124
125 if ( !gPythia->particleData.canDecay(decay_particle_pdg_code) ) {
126 LOG("Pythia8Decay", pWARN)
127 << " Particle of pdgid = " << decay_particle_pdg_code
128 << " can NOT be decayed by Pythia8";
129 return false;
130 }
131
132 // NOTE: Energy should be in GeV
133
134 gPythia->event.append( decay_particle_pdg_code, 1, 0, 0,
135 decay_particle_p4.Px(),
136 decay_particle_p4.Py(),
137 decay_particle_p4.Pz(),
138 decay_particle_p4.E(),
139 decay_particle_p4.M() );
140
141 // specify polarization, if any
142
143 // NOTE: while in Py8 polarization is a double variable ,
144 // in reality it's expected to be -1, 0., or 1 in case of "external" tau's,
145 // similar to LHA SPINUP; see Particle Decays, Hadron and Tau Decays in docs at
146 // https://pythia.org/manuals/pythia8305/Welcome.html
147 // so it's not able to handle anything like 0.99, thus we're rounding off
148 // gPythia->event.back().pol( round( std::cos( track.GetPolarization().angle( track.GetMomentumDirection() ) ) ) );
149 double a = decay_particle_p4.Vect().Angle(polz);
150 gPythia->event.back().pol( round( std::cos(a) ) );
151
152 int npart_before_decay = gPythia->event.size();
153
154 gPythia->next();
155
156 //gPythia->event.list();
157 //gPythia->stat();
158
159 int npart_after_decay = gPythia->event.size();
160
161 LOG("Pythia8Decay", pINFO) << "before " << npart_before_decay << " after " << npart_after_decay;
162
163 // the values of space coordinates from pythia are in mm.
164 // our conventions want it in fm
165 constexpr double space_scale = units::mm / units::fm ;
166
167 // the values of time coordinate from pythia is in mm/c.
168 // our conventions want it in ys
169 constexpr double time_scale = 1e21 * units::m / units::s ;
170
171
172 // create & fill up decay products
173 //
174 for ( int ip=npart_before_decay; ip<npart_after_decay; ++ip )
175 {
176
177 // only select final state decay products (direct or via subsequent decays);
178 // skip all others
179 //
180 // NOTE: in general, final state decays products will have
181 // positive status code between 91 and 99
182 // (in case such information could be of interest in the future)
183 //
184 if ( gPythia->event[ip].status() < 0 ) continue;
185
186 Pythia8::Event &fEvent = gPythia->event;
187
188 int daughter_pdg_code = fEvent[ip].id();
189
190 int mother1 = decay_particle_id;
191 int mother2 = -1;
192
193 int daughter1 = -1;
194 int daughter2 = -1;
195
196 bool is_hadron = pdg::IsHadron(daughter_pdg_code);
197 bool hadron_in_nuc = (in_nucleus && is_hadron && fRunBefHadroTransp);
198
199 GHepStatus_t daughter_status_code = (hadron_in_nuc) ?
201
202 GHepParticle mcp = GHepParticle(
203 daughter_pdg_code, // pdg
204 daughter_status_code, // status
205 mother1, // first parent
206 mother2, // second parent
207 daughter1, // first daughter
208 daughter2, // second daughter
209 fEvent[ip].px(), // px
210 fEvent[ip].py(), // py
211 fEvent[ip].pz(), // pz
212 fEvent[ip].e(), // e
213 fEvent[ip].xProd() * space_scale + decay_particle_x4.X(), // x
214 fEvent[ip].yProd() * space_scale + decay_particle_x4.Y(), // y
215 fEvent[ip].zProd() * space_scale + decay_particle_x4.Z(), // z
216 fEvent[ip].tProd() * time_scale + decay_particle_x4.T() // t
217 );
218
219 SLOG("Pythia8Decay", pINFO)
220 << "Adding daughter particle with PDG code = "
221 << daughter_pdg_code << ", m = " << mcp.Mass()
222 << " GeV, E = " << mcp.Energy() << " GeV)";
223
224 event->AddParticle(mcp);
225
226 }
227
228 // Update the event weight for each weighted particle decay
229 double weight = event->Weight() * fWeight;
230 event->SetWeight(weight);
231
232 // Mark input particle as a 'decayed state' & add its daughter links
233 decay_particle->SetStatus(kIStDecayedState);
234
235 return true;
236#else
237 LOG("Pythia8Decay", pFATAL)
238 << "calling GENIE/PYTHIA8 decay without enabling PYTHIA8";
239 gAbortingInErr = true;
240 std::exit(1);
241
242 return false;
243#endif
244
245}
#define pINFO
Definition Messenger.h:62
#define pFATAL
Definition Messenger.h:56
#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
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition Messenger.h:84
bool fRunBefHadroTransp
is invoked before or after FSI?
Definition Decayer.h:57
int Pdg(void) const
const TLorentzVector * P4(void) const
double Mass(void) const
Mass that corresponds to the PDG code.
const TLorentzVector * X4(void) const
void SetStatus(GHepStatus_t s)
double Energy(void) const
Get energy.
void GetPolarization(TVector3 &polz)
static Pythia8Singleton * Instance(void)
const double e
const double a
bool IsHadron(int pdgc)
Definition PDGUtils.cxx:392
static constexpr double m
Definition Units.h:71
static constexpr double mm
Definition Units.h:65
static constexpr double fm
Definition Units.h:75
static constexpr double s
Definition Units.h:95
@ kIStHadronInTheNucleus
Definition GHepStatus.h:37
@ kIStStableFinalState
Definition GHepStatus.h:30
@ kIStDecayedState
Definition GHepStatus.h:32
bool gAbortingInErr
Definition Messenger.cxx:34
enum genie::EGHepStatus GHepStatus_t

References a, e, genie::GHepParticle::Energy(), genie::units::fm, genie::Decayer::fRunBefHadroTransp, fWeight, genie::gAbortingInErr, genie::GHepParticle::GetPolarization(), genie::Pythia8Singleton::Instance(), genie::pdg::IsHadron(), genie::kIStDecayedState, genie::kIStHadronInTheNucleus, genie::kIStStableFinalState, LOG, genie::units::m, genie::GHepParticle::Mass(), genie::units::mm, genie::GHepParticle::P4(), genie::GHepParticle::Pdg(), pFATAL, pINFO, pWARN, genie::units::s, genie::GHepParticle::SetStatus(), SLOG, and genie::GHepParticle::X4().

Referenced by ProcessEventRecord().

◆ FindPythiaDecayChannel()

int Pythia8Decayer2023::FindPythiaDecayChannel ( int kc,
TDecayChannel * ch ) const
private

Definition at line 461 of file Pythia8Decayer2023.cxx.

462{
463 if(!dc) return -1;
464
465#ifdef __GENIE_PYTHIA8_ENABLED__
466
467 LOG("Pythia8Decay", pFATAL)
468 << "Pythia8Decayer2023::FindPythiaDecayChannel() not yet implemented";
469 gAbortingInErr = true;
470 std::exit(1);
471
472#else
473 LOG("Pythia8Decay", pFATAL)
474 << "calling GENIE/PYTHIA8 decay without enabling PYTHIA8";
475 gAbortingInErr = true;
476 std::exit(1);
477#endif
478
479 return -1;
480}

References genie::gAbortingInErr, LOG, and pFATAL.

◆ InhibitDecay()

void Pythia8Decayer2023::InhibitDecay ( int pdgc,
TDecayChannel * ch = 0 ) const
privatevirtual

Implements genie::Decayer.

Definition at line 295 of file Pythia8Decayer2023.cxx.

296{
297 if(! this->IsHandled(pdg_code)) return;
298
299#ifdef __GENIE_PYTHIA8_ENABLED__
300
301 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
302 bool known = gPythia->particleData.isParticle(pdg_code);
303 if ( ! known ) {
304 LOG("Pythia8Decay", pERROR)
305 << "Can not switch off decays of " << pdg_code
306 << " as it is UNKNOWN to Pythia8";
307 return;
308 }
309
310 auto pdentry = gPythia->particleData.particleDataEntryPtr(pdg_code);
311
312 int ifirst_chan = 0;
313 int ilast_chan = pdentry->sizeChannels() - 1;
314 if (ilast_chan < 0) {
315 LOG("Pythia8Decay", pDEBUG)
316 << "No available decay channels for " << pdg_code;
317 return;
318 }
319
320 if(dc) {
321 LOG("Pythia8Decay", pFATAL)
322 << "Switching off individual decay channels for particle = " << pdg_code
323 << " is NOT yet supported for Pythia8";
324 exit(42);
325 return;
326 }
327
328 //std::cout << "Before inhibiting channels " << ifirst_chan << "," << ilast_chan
329 // << " of " << pdg_code << std::endl;
330 //gPythia->particleData.list(pdg_code);
331
332 for (int ichan=ifirst_chan; ichan<=ilast_chan; ++ichan) {
333 int onMode = pdentry->channel(ichan).onMode();
334 // 4 possible modes: 0=off particle & antiparticle; 1=on both; 2=on particle only; 3=on antiparticle only
335 bool is_particle = (pdg_code>0);
336 // we're trying to turn off the mode
337 switch ( onMode ) {
338 case 0: // off for particles & antiparticles
339 break; // already off
340 case 1: // on for both particle & antiparticle
341 onMode = (is_particle ? 3 : 2);
342 break;
343 case 2: // on for particle only
344 onMode = (is_particle ? 0 : 2);
345 break;
346 case 3: // on for antiparticle only
347 onMode = (is_particle ? 3 : 0);
348 break;
349 }
350 pdentry->channel(ichan).onMode(onMode);
351 }
352
353 //std::cout << "After inhibiting channels " << ifirst_chan << "," << ilast_chan
354 // << " of " << pdg_code << std::endl;
355 //gPythia->particleData.list(pdg_code);
356
357#else
358 LOG("Pythia8Decay", pFATAL)
359 << "calling GENIE/PYTHIA8 decay without enabling PYTHIA8";
360 gAbortingInErr = true;
361 std::exit(1);
362#endif
363
364}
#define pERROR
Definition Messenger.h:59
#define pDEBUG
Definition Messenger.h:63

References genie::gAbortingInErr, genie::Pythia8Singleton::Instance(), IsHandled(), LOG, pDEBUG, pERROR, and pFATAL.

◆ Initialize()

void Pythia8Decayer2023::Initialize ( void ) const
private

Definition at line 247 of file Pythia8Decayer2023.cxx.

248{
249#ifdef __GENIE_PYTHIA8_ENABLED__
250 fWeight = 1.;
251
252 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
253 gPythia->readString("ProcessLevel:all = off");
254 gPythia->readString("ProcessLevel:resonanceDecays=on");
255
256 gPythia->readString("Print:quiet = on");
257
258 // sync GENIE and PYTHIA8 seeds
259 RandomGen * rnd = RandomGen::Instance();
260 long int seed = rnd->GetSeed();
261 gPythia->readString("Random::setSeed = on");
262 gPythia->settings.mode("Random:seed",seed);
263 LOG("Pythia8Decay", pINFO)
264 << "PYTHIA8 seed = " << gPythia->settings.mode("Random:seed");
265
266 gPythia->init();
267
268 // shut off decays of pi0's as we want Geant4 to handle them
269 // if other immediate decay products should be handled by Geant4,
270 // their respective decay modes should be shut off as well
271 //
272 gPythia->readString("111:onMode = off");
273
274#else
275 LOG("Pythia8Decay", pFATAL)
276 << "calling GENIE/PYTHIA8 decay without enabling PYTHIA8";
277 gAbortingInErr = true;
278 std::exit(1);
279#endif
280}
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
long int GetSeed(void) const
Definition RandomGen.h:82

References fWeight, genie::gAbortingInErr, genie::RandomGen::GetSeed(), genie::Pythia8Singleton::Instance(), genie::RandomGen::Instance(), LOG, pFATAL, and pINFO.

Referenced by Pythia8Decayer2023(), and Pythia8Decayer2023().

◆ IsHandled()

bool Pythia8Decayer2023::IsHandled ( int pdgc) const
privatevirtual

Implements genie::Decayer.

Definition at line 282 of file Pythia8Decayer2023.cxx.

283{
284// does not handle requests to decay baryon resonances
285
286 bool is_handled = (!utils::res::IsBaryonResonance(pdg_code));
287
288 LOG("Pythia8Decay", pDEBUG)
289 << "Can decay particle with PDG code = " << pdg_code
290 << "? " << ((is_handled)? "Yes" : "No");
291
292 return is_handled;
293}
bool IsBaryonResonance(int pdgc)
is input a baryon resonance?

References genie::utils::res::IsBaryonResonance(), LOG, and pDEBUG.

Referenced by InhibitDecay(), ProcessEventRecord(), and UnInhibitDecay().

◆ MatchDecayChannels()

bool Pythia8Decayer2023::MatchDecayChannels ( int ic,
TDecayChannel * ch ) const
private

Definition at line 482 of file Pythia8Decayer2023.cxx.

483{
484 // num. of daughters in the input TDecayChannel & the input PYTHIA ichannel
485 //int nd = dc->NDaughters();
486
487#ifdef __GENIE_PYTHIA8_ENABLED__
488
489 LOG("Pythia8Decay", pFATAL)
490 << "Pythia8Decayer2023::MatchDecayChannels() not yet implemented";
491 gAbortingInErr = true;
492 std::exit(1);
493
494#else
495 LOG("Pythia8Decay", pFATAL)
496 << "calling GENIE/PYTHIA8 decay without enabling PYTHIA8";
497 gAbortingInErr = true;
498 std::exit(1);
499#endif
500
501 return true;
502}

References genie::gAbortingInErr, LOG, and pFATAL.

◆ ProcessEventRecord()

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

Implements genie::EventRecordVisitorI.

Definition at line 54 of file Pythia8Decayer2023.cxx.

55{
56 LOG("ResonanceDecay", pINFO)
57 << "Running PYTHIA8 particle decayer "
58 << ((fRunBefHadroTransp) ? "*before*" : "*after*") << " FSI";
59
60 // Loop over particles, find unstable ones and decay them
61 TObjArrayIter piter(event);
62 GHepParticle * p = 0;
63 int ipos = -1;
64
65 while( (p = (GHepParticle *) piter.Next()) ) {
66 ipos++;
67
68 LOG("Pythia8Decay", pDEBUG) << "Checking: " << p->Name();
69
70 int pdg_code = p->Pdg();
71 GHepStatus_t status_code = p->Status();
72
73 if(!this->IsHandled (pdg_code)) continue;
74 if(!this->ToBeDecayed(pdg_code, status_code)) continue;
75
76 LOG("Pythia8Decay", pNOTICE)
77 << "Decaying unstable particle: " << p->Name();
78
79 bool decayed = this->Decay(ipos, event);
80 assert(decayed); // handle this more graciously and throw an exception
81 }
82
83 LOG("Pythia8Decay", pNOTICE)
84 << "Done finding & decaying unstable particles";
85}
#define pNOTICE
Definition Messenger.h:61
virtual bool ToBeDecayed(int pdgc, GHepStatus_t ist) const
Definition Decayer.cxx:51
string Name(void) const
Name that corresponds to the PDG code.
GHepStatus_t Status(void) const
bool Decay(int ip, GHepRecord *event) const

References Decay(), genie::Decayer::fRunBefHadroTransp, IsHandled(), LOG, genie::GHepParticle::Name(), pDEBUG, genie::GHepParticle::Pdg(), pINFO, pNOTICE, genie::GHepParticle::Status(), and genie::Decayer::ToBeDecayed().

◆ SumOfBranchingRatios()

double Pythia8Decayer2023::SumOfBranchingRatios ( int kc) const
private

Definition at line 436 of file Pythia8Decayer2023.cxx.

437{
438// Sum of branching ratios for enabled channels
439//
440 double sumbr=0.;
441
442#ifdef __GENIE_PYTHIA8_ENABLED__
443
444 LOG("Pythia8Decay", pFATAL)
445 << "Pythia8Decayer2023::SumOfBranchingRatios() not yet implemented";
446 gAbortingInErr = true;
447 std::exit(1);
448
449 LOG("Pythia8Decay", pINFO) << "Sum{BR} = " << sumbr;
450
451#else
452 LOG("Pythia8Decay", pFATAL)
453 << "calling GENIE/PYTHIA8 decay without enabling PYTHIA8";
454 gAbortingInErr = true;
455 std::exit(1);
456#endif
457
458 return sumbr;
459}

References genie::gAbortingInErr, LOG, pFATAL, and pINFO.

◆ UnInhibitDecay()

void Pythia8Decayer2023::UnInhibitDecay ( int pdgc,
TDecayChannel * ch = 0 ) const
privatevirtual

Implements genie::Decayer.

Definition at line 366 of file Pythia8Decayer2023.cxx.

367{
368 if(! this->IsHandled(pdg_code)) return;
369
370#ifdef __GENIE_PYTHIA8_ENABLED__
371
372 Pythia8::Pythia* gPythia = Pythia8Singleton::Instance()->Pythia8();
373 bool known = gPythia->particleData.isParticle(pdg_code);
374 if ( ! known ) {
375 LOG("Pythia8Decay", pERROR)
376 << "Can not switch on decays of " << pdg_code
377 << " as it is UNKNOWN to Pythia8";
378 return;
379 }
380
381 auto pdentry = gPythia->particleData.particleDataEntryPtr(pdg_code);
382
383 int ifirst_chan = 0;
384 int ilast_chan = pdentry->sizeChannels() - 1;
385 if (ilast_chan < 0) {
386 LOG("Pythia8Decay", pDEBUG)
387 << "No available decay channels for " << pdg_code;
388 return;
389 }
390
391 if(dc) {
392 LOG("Pythia8Decay", pFATAL)
393 << "Switching on individual decay channels for particle = " << pdg_code
394 << " is NOT yet supported for Pythia8";
395 exit(42);
396 return;
397 }
398
399 //std::cout << "Before uninhibiting channels " << ifirst_chan << "," << ilast_chan
400 // << " of " << pdg_code << std::endl;
401 //gPythia->particleData.list(pdg_code);
402
403 for (int ichan=ifirst_chan; ichan<=ilast_chan; ++ichan) {
404 int onMode = pdentry->channel(ichan).onMode();
405 // 4 possible modes: 0=off particle & antiparticle; 1=on both; 2=on particle only; 3=on antiparticle only
406 bool is_particle = (pdg_code>0);
407 // we're trying to turn off the mode
408 switch ( onMode ) {
409 case 0: // off for particles & antiparticles
410 onMode = (is_particle ? 2 : 3);
411 break;
412 case 1: // on for both particle & antiparticle
413 break; // already on
414 case 2: // on for particle only
415 onMode = (is_particle ? 2 : 0);
416 break;
417 case 3: // on for antiparticle only
418 onMode = (is_particle ? 0 : 3);
419 break;
420 }
421 pdentry->channel(ichan).onMode(onMode);
422 }
423
424 //std::cout << "After uninhibiting channels " << ifirst_chan << "," << ilast_chan
425 // << " of " << pdg_code << std::endl;
426 //gPythia->particleData.list(pdg_code);
427
428#else
429 LOG("Pythia8Decay", pFATAL)
430 << "calling GENIE/PYTHIA8 decay without enabling PYTHIA8";
431 gAbortingInErr = true;
432 std::exit(1);
433#endif
434}

References genie::gAbortingInErr, genie::Pythia8Singleton::Instance(), IsHandled(), LOG, pDEBUG, pERROR, and pFATAL.

Member Data Documentation

◆ fWeight

double genie::Pythia8Decayer2023::fWeight
mutableprivate

Definition at line 53 of file Pythia8Decayer2023.h.

Referenced by Decay(), and Initialize().


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