GENIEGenerator
Loading...
Searching...
No Matches
Pythia6Decayer2023.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//____________________________________________________________________________
10
11#include <vector>
12#include <cassert>
13
14#include <TClonesArray.h>
15#include <TLorentzVector.h>
16#include <TDecayChannel.h>
17#include <RVersion.h>
18
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>
23#else
24#include <TMCParticle6.h>
25#endif
26#endif
27
40
41using std::vector;
42
43using namespace genie;
44
45#ifdef __GENIE_PYTHIA6_ENABLED__
46// actual PYTHIA calls:
47extern "C" void py1ent_(int *, int *, double *, double *, double *);
48extern "C" void pydecy_(int *);
49#endif
50
51//____________________________________________________________________________
53Decayer("genie::Pythia6Decayer2023")
54{
55 this->Initialize();
56}
57//____________________________________________________________________________
59Decayer("genie::Pythia6Decayer2023", config)
60{
61 this->Initialize();
62}
63//____________________________________________________________________________
68//____________________________________________________________________________
70{
71 LOG("ResonanceDecay", pINFO)
72 << "Running PYTHIA6 particle decayer "
73 << ((fRunBefHadroTransp) ? "*before*" : "*after*") << " FSI";
74
75 // Loop over particles, find unstable ones and decay them
76 TObjArrayIter piter(event);
77 GHepParticle * p = 0;
78 int ipos = -1;
79
80 while( (p = (GHepParticle *) piter.Next()) ) {
81 ipos++;
82
83 LOG("Pythia6Decay", pDEBUG) << "Checking: " << p->Name();
84
85 int pdg_code = p->Pdg();
86 GHepStatus_t status_code = p->Status();
87
88 if(!this->IsHandled (pdg_code)) continue;
89 if(!this->ToBeDecayed(pdg_code, status_code)) continue;
90
91 LOG("Pythia6Decay", pNOTICE)
92 << "Decaying unstable particle: " << p->Name();
93
94 bool decayed = this->Decay(ipos, event);
95 assert(decayed); // handle this more graciously and throw an exception
96 }
97
98 LOG("Pythia6Decay", pNOTICE)
99 << "Done finding & decaying unstable particles";
100}
101//____________________________________________________________________________
102bool Pythia6Decayer2023::Decay(int decay_particle_id, GHepRecord * event) const
103{
104 fWeight = 1.; // reset previous decay weight
105
106#ifdef __GENIE_PYTHIA6_ENABLED__
107 // Get particle to be decayed
108 GHepParticle * decay_particle = event->Particle(decay_particle_id);
109 if(!decay_particle) return 0;
110
111 // Get the particle 4-momentum, 4-position and PDG code
112 TLorentzVector decay_particle_p4 = *(decay_particle->P4());
113 int decay_particle_pdg_code = decay_particle->Pdg();
114
115 // Convert to PYTHIA6 particle code and check whether decay is inhibited
116 int kc = fPythia->Pycomp(decay_particle_pdg_code);
117 int mdcy = fPythia->GetMDCY(kc, 1);
118 if(mdcy == 0) {
119 LOG("Pythia6Decay", pNOTICE)
120 << (PDGLibrary::Instance())->Find(decay_particle_pdg_code)->GetName()
121 << " decays are inhibited!";
122 return false;
123 }
124
125 // Get sub of BRs and compute weight if decay channels were inhibited
126 double sumbr = this->SumOfBranchingRatios(kc);
127 if(sumbr <= 0) {
128 LOG("Pythia6Decay", pWARN)
129 << "The sum of enabled "
130 << (PDGLibrary::Instance())->Find(decay_particle_pdg_code)->GetName()
131 << " decay channel branching ratios is non-positive!";
132 return false;
133 }
134 fWeight = 1./sumbr; // update weight to account for inhibited channels
135
136 // Run PYTHIA6 decay
137 int ip = 0;
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);
143
144 // Get decay products
145 fPythia->GetPrimaries();
146 TClonesArray * impl = (TClonesArray *) fPythia->ImportParticles("All");
147 if(!impl) {
148 LOG("Pythia6Decay", pWARN)
149 << "TPythia6::ImportParticles() returns a null list!";
150 return false;
151 }
152
153 // Copy the PYTHIA6 container to the GENIE event record
154
155 // Check whether the interaction is off a nuclear target or free nucleon
156 // Depending on whether this module is run before or after the hadron
157 // transport module it would affect the daughters status code
158 GHepParticle * target_nucleus = event->TargetNucleus();
159 bool in_nucleus = (target_nucleus!=0);
160
161 // the values of space coordinates from pythia are in mm.
162 // our conventions want it in fm
163 constexpr double space_scale = units::mm / units::fm ;
164
165 // the values of time coordinate from pythia is in mm/c.
166 // our conventions want it in ys
167 constexpr double time_scale = 1e21 * units::m / units::s ;
168
169 TMCParticle * p = 0;
170 TIter particle_iter(impl);
171 while( (p = (TMCParticle *) particle_iter.Next()) ) {
172 // Convert from TMCParticle to GHepParticle
174 p->GetKF(), // pdg
175 GHepStatus_t(p->GetKS()), // status
176 p->GetParent(), // first parent
177 0, // second parent
178 p->GetFirstChild(), // first daughter
179 p->GetLastChild(), // second daughter
180 p->GetPx(), // px
181 p->GetPy(), // py
182 p->GetPz(), // pz
183 p->GetEnergy(), // e
184 p->GetVx() * space_scale , // x
185 p->GetVy() * space_scale , // y
186 p->GetVz() * space_scale , // z
187 p->GetTime() * time_scale // t
188 );
189
190 if(mcp.Status()==kIStNucleonTarget) continue; // mother particle, already in GHEP
191
192 int daughter_pdg_code = mcp.Pdg();
193 SLOG("Pythia6Decay", pINFO)
194 << "Adding daughter particle wit PDG code = "
195 << daughter_pdg_code << ", m = " << mcp.Mass()
196 << " GeV, E = " << mcp.Energy() << " GeV)";
197
198 bool is_hadron = pdg::IsHadron(daughter_pdg_code);
199 bool hadron_in_nuc = (in_nucleus && is_hadron && fRunBefHadroTransp);
200
201 GHepStatus_t daughter_status_code = (hadron_in_nuc) ?
203
204 TLorentzVector daughter_p4(
205 mcp.Px(),mcp.Py(),mcp.Pz(),mcp.Energy());
206
207 event->AddParticle(
208 daughter_pdg_code, daughter_status_code,
209 decay_particle_id,-1,-1,-1,
210 daughter_p4, * mcp.X4() );
211 }
212
213 // Update the event weight for each weighted particle decay
214 double weight = event->Weight() * fWeight;
215 event->SetWeight(weight);
216
217 // Mark input particle as a 'decayed state' & add its daughter links
218 decay_particle->SetStatus(kIStDecayedState);
219
220 return true;
221#else
222 LOG("Pythia6Decay", pFATAL)
223 << "calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
224 gAbortingInErr = true;
225 std::exit(1);
226
227return false;
228#endif
229}
230//____________________________________________________________________________
232{
233#ifdef __GENIE_PYTHIA6_ENABLED__
234 fPythia = TPythia6::Instance();
235 fWeight = 1.;
236
237 // sync GENIE/PYTHIA6 seeds
239#else
240 LOG("Pythia6Decay", pFATAL)
241 << "calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
242 gAbortingInErr = true;
243 std::exit(1);
244#endif
245}
246//____________________________________________________________________________
247bool Pythia6Decayer2023::IsHandled(int pdg_code) const
248{
249// does not handle requests to decay baryon resonances
250
251 bool is_handled = (!utils::res::IsBaryonResonance(pdg_code));
252
253 LOG("Pythia6Decay", pDEBUG)
254 << "Can decay particle with PDG code = " << pdg_code
255 << "? " << ((is_handled)? "Yes" : "No");
256
257 return is_handled;
258}
259//____________________________________________________________________________
260void Pythia6Decayer2023::InhibitDecay(int pdg_code, TDecayChannel * dc) const
261{
262 if(! this->IsHandled(pdg_code)) return;
263
264#ifdef __GENIE_PYTHIA6_ENABLED__
265 int kc = fPythia->Pycomp(pdg_code);
266
267 if(!dc) {
268 LOG("Pythia6Decay", pINFO)
269 << "Switching OFF ALL decay channels for particle = " << pdg_code;
270 fPythia->SetMDCY(kc, 1,0);
271 return;
272 }
273
274 LOG("Pythia6Decay", pINFO)
275 << "Switching OFF decay channel = " << dc->Number()
276 << " for particle = " << pdg_code;
277
278 int ichannel = this->FindPythiaDecayChannel(kc, dc);
279 if(ichannel != -1) {
280 fPythia->SetMDME(ichannel,1,0); // switch-off
281 }
282#else
283 LOG("Pythia6Decay", pFATAL)
284 << "calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
285 gAbortingInErr = true;
286 std::exit(1);
287#endif
288
289}
290//____________________________________________________________________________
291void Pythia6Decayer2023::UnInhibitDecay(int pdg_code, TDecayChannel * dc) const
292{
293 if(! this->IsHandled(pdg_code)) return;
294
295#ifdef __GENIE_PYTHIA6_ENABLED__
296 int kc = fPythia->Pycomp(pdg_code);
297
298 if(!dc) {
299 LOG("Pythia6Decay", pINFO)
300 << "Switching ON all PYTHIA decay channels for particle = " << pdg_code;
301
302 fPythia->SetMDCY(kc, 1,1);
303
304 int first_channel = fPythia->GetMDCY(kc,2);
305 int last_channel = fPythia->GetMDCY(kc,2) + fPythia->GetMDCY(kc,3) - 1;
306
307 for(int ichannel = first_channel;
308 ichannel < last_channel; ichannel++) {
309 fPythia->SetMDME(ichannel,1,1); // switch-on
310 }
311 return;
312 }//!dc
313
314 LOG("Pythia6Decay", pINFO)
315 << "Switching OFF decay channel = " << dc->Number()
316 << " for particle = " << pdg_code;
317
318 int ichannel = this->FindPythiaDecayChannel(kc, dc);
319 if(ichannel != -1) {
320 fPythia->SetMDME(ichannel,1,1); // switch-on
321 }
322#else
323 LOG("Pythia6Decay", pFATAL)
324 << "calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
325 gAbortingInErr = true;
326 std::exit(1);
327#endif
328}
329//____________________________________________________________________________
331{
332// Sum of branching ratios for enabled channels
333//
334 double sumbr=0.;
335
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;
339
340 bool has_inhibited_channels=false;
341
342 // loop over pythia decay channels
343 for(int ichannel = first_channel;
344 ichannel < last_channel; ichannel++) {
345
346 bool enabled = (fPythia->GetMDME(ichannel,1) == 1);
347 if (!enabled) {
348 has_inhibited_channels = true;
349 } else {
350 sumbr += fPythia->GetBRAT(ichannel);
351 }
352 }
353
354 if(!has_inhibited_channels) return 1.;
355 LOG("Pythia6Decay", pINFO) << "Sum{BR} = " << sumbr;
356#else
357 LOG("Pythia6Decay", pFATAL)
358 << "calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
359 gAbortingInErr = true;
360 std::exit(1);
361#endif
362
363 return sumbr;
364}
365//____________________________________________________________________________
366int Pythia6Decayer2023::FindPythiaDecayChannel(int kc, TDecayChannel* dc) const
367{
368 if(!dc) return -1;
369
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;
373
374 bool found_match = false;
375
376 // loop over pythia decay channels
377 for(int ichannel = first_channel;
378 ichannel < last_channel; ichannel++) {
379
380 // does the current pythia channel matches the input TDecayChannel?
381 LOG("Pythia6Decay", pINFO)
382 << "\nComparing PYTHIA's channel = " << ichannel
383 << " with TDecayChannel = " << dc->Number();
384
385 found_match = this->MatchDecayChannels(ichannel,dc);
386 if(found_match) {
387 LOG("Pythia6Decay", pNOTICE)
388 << " ** TDecayChannel id = " << dc->Number()
389 << " corresponds to PYTHIA6 channel id = " << ichannel;
390 return ichannel;
391 }//match?
392 }//loop pythia decay ch.
393
394 LOG("Pythia6Decay", pWARN)
395 << " ** No PYTHIA6 decay channel match found for "
396 << "TDecayChannel id = " << dc->Number();
397#else
398 LOG("Pythia6Decay", pFATAL)
399 << "calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
400 gAbortingInErr = true;
401 std::exit(1);
402#endif
403
404
405 return -1;
406}
407//____________________________________________________________________________
408bool Pythia6Decayer2023::MatchDecayChannels(int ichannel, TDecayChannel* dc) const
409{
410 // num. of daughters in the input TDecayChannel & the input PYTHIA ichannel
411 int nd = dc->NDaughters();
412
413#ifdef __GENIE_PYTHIA6_ENABLED__
414 int py_nd = 0;
415 for (int i = 1; i <= 5; i++) {
416 if(fPythia->GetKFDP(ichannel,i) != 0) py_nd++;
417 }
418
419 LOG("Pythia6Decay", pDEBUG)
420 << "NDaughters: PYTHIA = " << py_nd << ", ROOT's TDecayChannel = " << nd;
421
422 if(nd != py_nd) return false;
423
424 //
425 // if the two channels have the same num. of daughters, then compare them
426 //
427
428 // store decay daughters for the input TDecayChannel
429 vector<int> dc_daughter(nd);
430 int k=0;
431 for( ; k < nd; k++) {
432 dc_daughter[k] = dc->DaughterPdgCode(k);
433 }
434 // store decay daughters for the input PYTHIA's ichannel
435 vector<int> py_daughter(nd);
436 k=0;
437 for(int i = 1; i <= 5; i++) {
438 if(fPythia->GetKFDP(ichannel,i) == 0) continue;
439 py_daughter[k] = fPythia->GetKFDP(ichannel,i);
440 k++;
441 }
442
443 // sort both daughter lists
444 sort( dc_daughter.begin(), dc_daughter.end() );
445 sort( py_daughter.begin(), py_daughter.end() );
446
447 // compare
448 for(int i = 0; i < nd; i++) {
449 if(dc_daughter[i] != py_daughter[i]) return false;
450 }
451#else
452 LOG("Pythia6Decay", pFATAL)
453 << "calling GENIE/PYTHIA6 decay without enabling PYTHIA6";
454 gAbortingInErr = true;
455 std::exit(1);
456#endif
457
458 return true;
459}
460//____________________________________________________________________________
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#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
#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
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?
Definition Decayer.h:57
virtual bool ToBeDecayed(int pdgc, GHepStatus_t ist) const
Definition Decayer.cxx:51
STDHEP-like event record entry that can fit a particle or a nucleus.
string Name(void) const
Name that corresponds to the PDG code.
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 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.
Definition GHepRecord.h:45
static PDGLibrary * Instance(void)
double SumOfBranchingRatios(int kc) const
void UnInhibitDecay(int pdgc, TDecayChannel *ch=0) const
bool Decay(int ip, GHepRecord *event) const
int FindPythiaDecayChannel(int kc, TDecayChannel *ch) const
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.
Definition RandomGen.cxx:74
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
bool IsBaryonResonance(int pdgc)
is input a baryon resonance?
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
@ kIStHadronInTheNucleus
Definition GHepStatus.h:37
@ kIStStableFinalState
Definition GHepStatus.h:30
@ kIStDecayedState
Definition GHepStatus.h:32
@ kIStNucleonTarget
Definition GHepStatus.h:34
bool gAbortingInErr
Definition Messenger.cxx:34
enum genie::EGHepStatus GHepStatus_t