14#include <TClonesArray.h>
15#include <TLorentzVector.h>
16#include <TDecayChannel.h>
19#include "Framework/Conventions/GBuild.h"
20#ifdef __GENIE_PYTHIA6_ENABLED__
21#if ROOT_VERSION_CODE >= ROOT_VERSION(5,15,6)
22#include <TMCParticle.h>
24#include <TMCParticle6.h>
45#ifdef __GENIE_PYTHIA6_ENABLED__
47extern "C" void py1ent_(
int *,
int *,
double *,
double *,
double *);
48extern "C" void pydecy_(
int *);
53Decayer(
"genie::Pythia6Decayer2023")
59Decayer(
"genie::Pythia6Decayer2023", config)
72 <<
"Running PYTHIA6 particle decayer "
76 TObjArrayIter piter(event);
85 int pdg_code = p->
Pdg();
89 if(!this->
ToBeDecayed(pdg_code, status_code))
continue;
92 <<
"Decaying unstable particle: " << p->
Name();
94 bool decayed = this->
Decay(ipos, event);
99 <<
"Done finding & decaying unstable particles";
106#ifdef __GENIE_PYTHIA6_ENABLED__
108 GHepParticle * decay_particle =
event->Particle(decay_particle_id);
109 if(!decay_particle)
return 0;
112 TLorentzVector decay_particle_p4 = *(decay_particle->
P4());
113 int decay_particle_pdg_code = decay_particle->
Pdg();
116 int kc = fPythia->Pycomp(decay_particle_pdg_code);
117 int mdcy = fPythia->GetMDCY(kc, 1);
121 <<
" decays are inhibited!";
129 <<
"The sum of enabled "
131 <<
" decay channel branching ratios is non-positive!";
138 double E = decay_particle_p4.Energy();
139 double theta = decay_particle_p4.Theta();
140 double phi = decay_particle_p4.Phi();
141 fPythia->SetMSTJ(22,1);
142 py1ent_(&ip, &decay_particle_pdg_code, &E, &theta, &phi);
145 fPythia->GetPrimaries();
146 TClonesArray * impl = (TClonesArray *) fPythia->ImportParticles(
"All");
149 <<
"TPythia6::ImportParticles() returns a null list!";
159 bool in_nucleus = (target_nucleus!=0);
170 TIter particle_iter(impl);
171 while( (p = (TMCParticle *) particle_iter.Next()) ) {
184 p->GetVx() * space_scale ,
185 p->GetVy() * space_scale ,
186 p->GetVz() * space_scale ,
187 p->GetTime() * time_scale
192 int daughter_pdg_code = mcp.
Pdg();
194 <<
"Adding daughter particle wit PDG code = "
195 << daughter_pdg_code <<
", m = " << mcp.
Mass()
196 <<
" GeV, E = " << mcp.
Energy() <<
" GeV)";
204 TLorentzVector daughter_p4(
208 daughter_pdg_code, daughter_status_code,
209 decay_particle_id,-1,-1,-1,
210 daughter_p4, * mcp.
X4() );
214 double weight =
event->Weight() *
fWeight;
215 event->SetWeight(weight);
223 <<
"calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
233#ifdef __GENIE_PYTHIA6_ENABLED__
234 fPythia = TPythia6::Instance();
241 <<
"calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
254 <<
"Can decay particle with PDG code = " << pdg_code
255 <<
"? " << ((is_handled)?
"Yes" :
"No");
264#ifdef __GENIE_PYTHIA6_ENABLED__
265 int kc = fPythia->Pycomp(pdg_code);
269 <<
"Switching OFF ALL decay channels for particle = " << pdg_code;
270 fPythia->SetMDCY(kc, 1,0);
275 <<
"Switching OFF decay channel = " << dc->Number()
276 <<
" for particle = " << pdg_code;
280 fPythia->SetMDME(ichannel,1,0);
284 <<
"calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
295#ifdef __GENIE_PYTHIA6_ENABLED__
296 int kc = fPythia->Pycomp(pdg_code);
300 <<
"Switching ON all PYTHIA decay channels for particle = " << pdg_code;
302 fPythia->SetMDCY(kc, 1,1);
304 int first_channel = fPythia->GetMDCY(kc,2);
305 int last_channel = fPythia->GetMDCY(kc,2) + fPythia->GetMDCY(kc,3) - 1;
307 for(
int ichannel = first_channel;
308 ichannel < last_channel; ichannel++) {
309 fPythia->SetMDME(ichannel,1,1);
315 <<
"Switching OFF decay channel = " << dc->Number()
316 <<
" for particle = " << pdg_code;
320 fPythia->SetMDME(ichannel,1,1);
324 <<
"calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
336#ifdef __GENIE_PYTHIA6_ENABLED__
337 int first_channel = fPythia->GetMDCY(kc,2);
338 int last_channel = fPythia->GetMDCY(kc,2) + fPythia->GetMDCY(kc,3) - 1;
340 bool has_inhibited_channels=
false;
343 for(
int ichannel = first_channel;
344 ichannel < last_channel; ichannel++) {
346 bool enabled = (fPythia->GetMDME(ichannel,1) == 1);
348 has_inhibited_channels =
true;
350 sumbr += fPythia->GetBRAT(ichannel);
354 if(!has_inhibited_channels)
return 1.;
355 LOG(
"Pythia6Decay",
pINFO) <<
"Sum{BR} = " << sumbr;
358 <<
"calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
370#ifdef __GENIE_PYTHIA6_ENABLED__
371 int first_channel = fPythia->GetMDCY(kc,2);
372 int last_channel = fPythia->GetMDCY(kc,2) + fPythia->GetMDCY(kc,3) - 1;
374 bool found_match =
false;
377 for(
int ichannel = first_channel;
378 ichannel < last_channel; ichannel++) {
382 <<
"\nComparing PYTHIA's channel = " << ichannel
383 <<
" with TDecayChannel = " << dc->Number();
388 <<
" ** TDecayChannel id = " << dc->Number()
389 <<
" corresponds to PYTHIA6 channel id = " << ichannel;
395 <<
" ** No PYTHIA6 decay channel match found for "
396 <<
"TDecayChannel id = " << dc->Number();
399 <<
"calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
411 int nd = dc->NDaughters();
413#ifdef __GENIE_PYTHIA6_ENABLED__
415 for (
int i = 1; i <= 5; i++) {
416 if(fPythia->GetKFDP(ichannel,i) != 0) py_nd++;
420 <<
"NDaughters: PYTHIA = " << py_nd <<
", ROOT's TDecayChannel = " << nd;
422 if(nd != py_nd)
return false;
429 vector<int> dc_daughter(nd);
431 for( ; k < nd; k++) {
432 dc_daughter[k] = dc->DaughterPdgCode(k);
435 vector<int> py_daughter(nd);
437 for(
int i = 1; i <= 5; i++) {
438 if(fPythia->GetKFDP(ichannel,i) == 0)
continue;
439 py_daughter[k] = fPythia->GetKFDP(ichannel,i);
444 sort( dc_daughter.begin(), dc_daughter.end() );
445 sort( py_daughter.begin(), py_daughter.end() );
448 for(
int i = 0; i < nd; i++) {
449 if(dc_daughter[i] != py_daughter[i])
return false;
453 <<
"calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
bool fRunBefHadroTransp
is invoked before or after FSI?
virtual bool ToBeDecayed(int pdgc, GHepStatus_t ist) const
STDHEP-like event record entry that can fit a particle or a nucleus.
string Name(void) const
Name that corresponds to the PDG code.
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 Px(void) const
Get Px.
double Pz(void) const
Get Pz.
double Py(void) const
Get Py.
double Energy(void) const
Get energy.
GHepStatus_t Status(void) const
GENIE's GHEP MC event record.
static PDGLibrary * Instance(void)
void Initialize(void) const
double SumOfBranchingRatios(int kc) const
void UnInhibitDecay(int pdgc, TDecayChannel *ch=0) const
bool Decay(int ip, GHepRecord *event) const
bool IsHandled(int pdgc) const
int FindPythiaDecayChannel(int kc, TDecayChannel *ch) const
virtual ~Pythia6Decayer2023()
void InhibitDecay(int pdgc, TDecayChannel *ch=0) const
bool MatchDecayChannels(int ic, TDecayChannel *ch) const
void ProcessEventRecord(GHepRecord *event) const
static RandomGen * Instance()
Access instance.
static constexpr double m
static constexpr double mm
static constexpr double fm
static constexpr double s
bool IsBaryonResonance(int pdgc)
is input a baryon resonance?
THE MAIN GENIE PROJECT NAMESPACE
enum genie::EGHepStatus GHepStatus_t