GENIEGenerator
Loading...
Searching...
No Matches
GLRESWdecPythia8.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 Alfonso Garcia <aagarciasoto \at km3net.de>
7 IFIC (Valencia)
8*/
9//____________________________________________________________________________
10
12
13//___________________________________________________________________________
15EventRecordVisitorI("genie::GLRESWdecPythia8")
16{
17 this->Initialize();
18 born = new Born();
19}
20//___________________________________________________________________________
22EventRecordVisitorI("genie::GLRESWdecPythia8", config)
23{
24 this->Initialize();
25 born = new Born();
26}
27//___________________________________________________________________________
32//____________________________________________________________________________
34{
35
36 if(!this->Wdecay(event)) {
37 LOG("GLRESGenerator", pWARN) << "W decayer failed!";
38 event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
40 exception.SetReason("Could not simulate the W decay system");
41 exception.SwitchOnFastForward();
42 throw exception;
43 return;
44 }
45
46}
47//___________________________________________________________________________
49#ifdef __GENIE_PYTHIA8_ENABLED__
50 event // avoid unused variable warning if PYTHIA8 is not enabled
51#endif
52) const
53{
54
55#ifdef __GENIE_PYTHIA8_ENABLED__
56
57 Interaction * interaction = event->Summary();
58 const InitialState & init_state = interaction->InitState();
59
60 //incoming v & struck particle & remnant nucleus
61 GHepParticle * nu = event->Probe();
62 GHepParticle * el = event->HitElectron();
63
64 GHepParticle * target = event -> TargetNucleus();
65 if(target) event->AddParticle(target->Pdg(), kIStFinalStateNuclearRemnant, 1,-1,-1,-1, *(target->P4()), *(target->X4()) );
66
67 TVector3 unit_nu = nu->P4()->Vect().Unit();
68
69 long double mlout = interaction->FSPrimLepton()->Mass(); //mass of charged lepton
70 long double mlin = kElectronMass; //mass of incoming charged lepton
71
72 long double Enuin = init_state.ProbeE(kRfLab);
73 long double s = born->GetS(mlin,Enuin);
74
75 long double n1 = interaction->Kine().GetKV(kKVn1);
76 long double n2 = interaction->Kine().GetKV(kKVn2);
77
78 long double costhCM = n1;
79 long double sinthCM = sqrtl(1-costhCM*costhCM);
80
81 long double t = born->GetT(mlin,mlout,s,n1);
82 long double zeta = born->GetReAlpha()/kPi*(2.0*logl(sqrtl(-t)/kElectronMass)-1.0);
83 long double omx = powl(n2, 1.0/zeta );
84 long double s_r = s*( 1.-omx );
85
86 // Boost velocity CM -> LAB
87 long double EnuinCM = (s_r-mlin*mlin)/sqrtl(s_r)/2.;
88 long double beta = (powl(Enuin,2)-powl(EnuinCM,2))/(powl(Enuin,2)+powl(EnuinCM,2));
89
90 // Final state primary lepton PDG code
91 int pdgl = interaction->FSPrimLeptonPdg();
92 assert(pdgl!=0);
93
94 if ( pdg::IsElectron(TMath::Abs(pdgl)) || pdg::IsMuon(TMath::Abs(pdgl)) || pdg::IsTau(TMath::Abs(pdgl)) ) {
95
96 long double ElpoutCM = (s_r+mlout*mlout)/sqrtl(s_r)/2.;
97 long double EnuoutCM = (s_r-mlout*mlout)/sqrtl(s_r)/2.;
98 LongLorentzVector p4_lpout( 0., EnuoutCM*sinthCM, EnuoutCM*costhCM, ElpoutCM );
99 LongLorentzVector p4_nuout( 0., -EnuoutCM*sinthCM, -EnuoutCM*costhCM, EnuoutCM );
100
101 p4_lpout.BoostZ(beta);
102 p4_nuout.BoostZ(beta);
103
104 TLorentzVector p4lp_o( (double)p4_lpout.Px(), (double)p4_lpout.Py(), (double)p4_lpout.Pz(), (double)p4_lpout.E() );
105 TLorentzVector p4nu_o( (double)p4_nuout.Px(), (double)p4_nuout.Py(), (double)p4_nuout.Pz(), (double)p4_nuout.E() );
106
107 // Randomize transverse components
109 double phi = 2* kPi * rnd->RndLep().Rndm();
110 p4lp_o.RotateZ(phi);
111 p4nu_o.RotateZ(phi);
112
113 //rotate from LAB=[0,0,Ev,Ev]->[px,py,pz,E]
114 p4lp_o.RotateUz(unit_nu);
115 p4nu_o.RotateUz(unit_nu);
116
117 int pdgvout = 0;
118 if ( pdg::IsElectron(pdgl) ) pdgvout = kPdgAntiNuE;
119 else if ( pdg::IsPositron(pdgl) ) pdgvout = kPdgNuE;
120 else if ( pdg::IsMuon(pdgl) ) pdgvout = kPdgAntiNuMu;
121 else if ( pdg::IsAntiMuon(pdgl) ) pdgvout = kPdgNuMu;
122 else if ( pdg::IsTau(pdgl) ) pdgvout = kPdgAntiNuTau;
123 else if ( pdg::IsAntiTau(pdgl) ) pdgvout = kPdgNuTau;
124
125 int pdgboson = pdg::IsNeutrino(init_state.ProbePdg()) ? kPdgWP : kPdgWM;
126
127 // Create a GHepParticle and add it to the event record
128 event->AddParticle( pdgboson, kIStDecayedState, 0, -1, 5, 6, *nu->P4()+*el->P4(), *(nu->X4()) ); //W [mothers: nuebar_in,e_in][daugthers: nulbar_out,lp_out]
129 event->AddParticle( pdgl, kIStStableFinalState, 4, -1, -1, -1, p4lp_o, *(nu->X4()) );
130 event->AddParticle( pdgvout, kIStStableFinalState, 4, -1, -1, -1, p4nu_o, *(nu->X4()) );
131 event->Summary()->KinePtr()->SetFSLeptonP4(p4lp_o);
132
133 }
134 else {
135
136 //We have to initialize Pythia8 everytime we generate a new event
137 //And we have to change the seed because it would regenerate the same event
138 //if we dont do it.
140 fPythia->settings.mode("Random:seed", rnd->RndLep().Integer(100000000));
141 fPythia->settings.parm("Beams:eCM",sqrtl(s_r));
142 fPythia->init();
143 fPythia->next();
144
145 // fPythia->event.list();
146 // fPythia->stat();
147
148 Pythia8::Event &fEvent = fPythia->event;
149 int np = fEvent.size();
150 assert(np>0);
151
152 for (int i = 5; i < np; i++) { // ignore first entry -> system + input particles
153
154 int pdgc = fEvent[i].id();
155 if (!PDGLibrary::Instance()->Find(pdgc)) continue; // some intermediate particles not part of genie tables
156
157 int ks = fEvent[i].status();
158
159 LongLorentzVector p4longo(fEvent[i].px(), fEvent[i].py(), fEvent[i].pz(), fEvent[i].e());
160 p4longo.BoostZ(beta);
161
162 TLorentzVector p4o( (double)p4longo.Px(), (double)p4longo.Py(), (double)p4longo.Pz(), (double)p4longo.E() );
163 p4o.RotateUz(unit_nu);
164
165 double massPDG = PDGLibrary::Instance()->Find(pdgc)->Mass();
166 if ( ks>0 && p4o.E()<massPDG ) {
167 LOG("GLRESGenerator", pWARN) << "Putting at rest one stable particle generated by PYTHIA because E < m";
168 LOG("GLRESGenerator", pWARN) << "PDG = " << pdgc << " // State = " << ks;
169 LOG("GLRESGenerator", pWARN) << "E = " << p4o.E() << " // |p| = " << TMath::Sqrt(p4o.P());
170 LOG("GLRESGenerator", pWARN) << "p = [ " << p4o.Px() << " , " << p4o.Py() << " , " << p4o.Pz() << " ]";
171 LOG("GLRESGenerator", pWARN) << "m = " << p4o.M() << " // mpdg = " << massPDG;
172 p4o.SetXYZT(0,0,0,massPDG);
173 }
174
175 // copy final state particles to the event record
177
178 // fix numbering scheme used for mother/daughter assignments
179 int firstmother = -1;
180 int lastmother = -1;
181 int firstchild = -1;
182 int lastchild = -1;
183
184 if ( fEvent[i].mother1() == 3 ) { //produced W boson: mother will be the incoming neutrino
185 firstmother = 0;
186 firstchild = fEvent[i].daughter1() - 1;
187 lastchild = fEvent[i].daughter2() - 1;
188 }
189 else { //rest
190 firstmother = fEvent[i].mother1() - 1; //shift to match boson position
191 firstchild = fEvent[i].daughter1() - 1;
192 lastchild = fEvent[i].daughter2() - 1;
193 }
194
195 double vx = nu->X4()->X() + fEvent[i].xProd()*1e12; //pythia gives position in [mm] while genie uses [fm]
196 double vy = nu->X4()->Y() + fEvent[i].yProd()*1e12;
197 double vz = nu->X4()->Z() + fEvent[i].zProd()*1e12;
198 double vt = nu->X4()->T() + fEvent[i].tProd()*(units::millimeter/units::second);
199 TLorentzVector pos( vx, vy, vz, vt );
200
201 event->AddParticle(pdgc, ist, firstmother, lastmother, firstchild, lastchild, p4o, pos );
202 }
203
204 }
205
206 return true;
207#else
208 return false;
209#endif
210
211}
212//___________________________________________________________________________
218//____________________________________________________________________________
224//____________________________________________________________________________
226{
227
228#ifdef __GENIE_PYTHIA8_ENABLED__
229 GetParam( "SSBarSuppression", fSSBarSuppression );
230 GetParam( "GaussianPt2", fGaussianPt2 );
231 GetParam( "NonGaussianPt2Tail", fNonGaussianPt2Tail );
232 GetParam( "RemainingEnergyCutoff", fRemainingECutoff );
233 GetParam( "DiQuarkSuppression", fDiQuarkSuppression );
234 GetParam( "LightVMesonSuppression", fLightVMesonSuppression );
235 GetParam( "SVMesonSuppression", fSVMesonSuppression );
236 GetParam( "Lunda", fLunda );
237 GetParam( "Lundb", fLundb );
238 GetParam( "LundaDiq", fLundaDiq );
239
240 fPythia->settings.parm("StringFlav:probStoUD", fSSBarSuppression);
241 fPythia->settings.parm("Diffraction:primKTwidth", fGaussianPt2);
242 fPythia->settings.parm("StringPT:enhancedFraction", fNonGaussianPt2Tail);
243 fPythia->settings.parm("StringFragmentation:stopMass", fRemainingECutoff);
244 fPythia->settings.parm("StringFlav:probQQtoQ", fDiQuarkSuppression);
245 fPythia->settings.parm("StringFlav:mesonUDvector", fLightVMesonSuppression);
246 fPythia->settings.parm("StringFlav:mesonSvector", fSVMesonSuppression);
247 fPythia->settings.parm("StringZ:aLund", fLunda);
248 fPythia->settings.parm("StringZ:bLund", fLundb);
249 fPythia->settings.parm("StringZ:aExtraDiquark", fLundaDiq);
250
251 // Same default mass of the W boson in pythia8 and genie, so no need to change in pythia8
252 // No problem with energy conservation W and top decays, so no need to set the width to 0
253
254 // Important for not decaying too early!!!
255 // Copy of the method from Decayer class to have a consistent treatment.
256 // The idea is to not decay in this stage what we dont want to decay in Decayer.
257 // Hadronization -> Decayer
258 //------------------------------
259 // Decay -> Decay
260 // Undecay -> Undecay
261 // Undecay -> Decay
262 // Decay -> Undecay -> This is the problematic one
263 // So we loop over all the particles for which we inhibit decay in the config file
264 // for the Decayer stage and we inhibit it also at hadronization.
265 // This is not need in LeptonHadronization and PhotonCOH because they use the
266 // Pythia8Singleton class which already defines the decays consistently
267 //
268 RgKeyList klist = GetConfig().FindKeys("DecayParticleWithCode=");
269 RgKeyList::const_iterator kiter = klist.begin();
270 for( ; kiter != klist.end(); ++kiter) {
271 RgKey key = *kiter;
272 bool decay = GetConfig().GetBool(key);
273 vector<string> kv = utils::str::Split(key,"=");
274 int pdg_code = atoi(kv[1].c_str());
275 if(!decay) {
276 LOG("GLRESGenerator", pDEBUG) << "Configured to inhibit decays for " << pdg_code;
277 auto pdentry = fPythia->particleData.particleDataEntryPtr(pdg_code);
278 for (int ichan=0; ichan<=pdentry->sizeChannels()-1; ++ichan) {
279 int onMode = pdentry->channel(ichan).onMode();
280 bool is_particle = (pdg_code>0);
281 switch ( onMode ) {
282 case 0: // off for particles & antiparticles
283 break; // already off
284 case 1: // on for both particle & antiparticle
285 onMode = (is_particle ? 3 : 2);
286 break;
287 case 2: // on for particle only
288 onMode = (is_particle ? 0 : 2);
289 break;
290 case 3: // on for antiparticle only
291 onMode = (is_particle ? 3 : 0);
292 break;
293 }
294 pdentry->channel(ichan).onMode(onMode);
295 }
296 }
297 }
298
299 //Important: We dont initialize Pythia8 here because it must be done event by event. See above
300#endif
301
302}
303//____________________________________________________________________________
305{
306
307#ifdef __GENIE_PYTHIA8_ENABLED__
308 fPythia = new Pythia8::Pythia();
309
310 fPythia->readString("Print:quiet = on");
311 fPythia->readString("Random:setSeed = on");
312
313 //One cool feature of having independent Pythia8 instances
314 //We can define our intial state with the proper setting (i.e. disabling decays)
315 //without affecting other classes of GENIE that also use Pythia8
316 fPythia->readString("WeakSingleBoson:ffbar2ffbar(s:W) = on");
317 fPythia->readString("PDF:lepton = off");
318 fPythia->readString("24:onMode = off");
319 fPythia->readString("24:onIfAny = 1 2 3 4 5"); //enable W->hadron only
320 fPythia->readString("Beams:idA = -12");
321 fPythia->readString("Beams:idB = 11");
322#endif
323
324}
#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
string RgKey
virtual const Registry & GetConfig(void) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
virtual void Configure(const Registry &config)
Definition Algorithm.cxx:62
Born level nu-electron cross section.
Definition Born.h:26
STDHEP-like event record entry that can fit a particle or a nucleus.
int Pdg(void) const
const TLorentzVector * P4(void) const
const TLorentzVector * X4(void) const
GENIE's GHEP MC event record.
Definition GHepRecord.h:45
double fDiQuarkSuppression
di-quark suppression parameter
double fLundaDiq
adjustment of Lund a for di-quark
double fLightVMesonSuppression
light vector meson suppression
double fRemainingECutoff
remaining E cutoff stopping fragmentation
double fNonGaussianPt2Tail
non gaussian pt2 tail parameterization
void Configure(const Registry &config)
double fLundb
Lund b parameter.
double fLunda
Lund a parameter.
double fSVMesonSuppression
strange vector meson suppression
void ProcessEventRecord(GHepRecord *event) const
double fGaussianPt2
gaussian pt2 distribution width
double fSSBarSuppression
ssbar suppression
bool Wdecay(GHepRecord *event) const
Initial State information.
int ProbePdg(void) const
double ProbeE(RefFrame_t rf) const
Summary information for an interaction.
Definition Interaction.h:56
const Kinematics & Kine(void) const
Definition Interaction.h:71
int FSPrimLeptonPdg(void) const
final state primary lepton pdg
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
const InitialState & InitState(void) const
Definition Interaction.h:69
double GetKV(KineVar_t kv) const
static PDGLibrary * Instance(void)
TParticlePDG * Find(int pdgc, bool must_exist=true)
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
TRandom3 & RndLep(void) const
rnd number generator used by final state primary lepton generators
Definition RandomGen.h:62
A registry. Provides the container for algorithm configuration parameters.
Definition Registry.h:65
RgKeyList FindKeys(RgKey key_part) const
create list with all keys containing 'key_part'
Definition Registry.cxx:840
RgBool GetBool(RgKey key) const
Definition Registry.cxx:460
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
const double e
bool IsTau(int pdgc)
Definition PDGUtils.cxx:208
bool IsNeutrino(int pdgc)
Definition PDGUtils.cxx:110
bool IsPositron(int pdgc)
Definition PDGUtils.cxx:193
bool IsElectron(int pdgc)
Definition PDGUtils.cxx:188
bool IsAntiTau(int pdgc)
Definition PDGUtils.cxx:213
bool IsAntiMuon(int pdgc)
Definition PDGUtils.cxx:203
bool IsMuon(int pdgc)
Definition PDGUtils.cxx:198
static constexpr double millimeter
Definition Units.h:41
static constexpr double second
Definition Units.h:37
Simple functions for loading and reading nucleus dependent keys from config files.
vector< string > Split(string input, string delim)
@ kIStFinalStateNuclearRemnant
Definition GHepStatus.h:38
@ kIStDISPreFragmHadronicState
Definition GHepStatus.h:35
@ kIStStableFinalState
Definition GHepStatus.h:30
@ kIStDecayedState
Definition GHepStatus.h:32
vector< RgKey > RgKeyList
Definition Registry.h:50
const int kPdgWM
Definition PDGCodes.h:192
const int kPdgAntiNuE
Definition PDGCodes.h:29
enum genie::EGHepStatus GHepStatus_t
const int kPdgAntiNuTau
Definition PDGCodes.h:33
@ kKVn2
Definition KineVar.h:62
@ kKVn1
Definition KineVar.h:61
const int kPdgNuE
Definition PDGCodes.h:28
const int kPdgNuTau
Definition PDGCodes.h:32
@ kRfLab
Definition RefFrame.h:26
const int kPdgAntiNuMu
Definition PDGCodes.h:31
const int kPdgWP
Definition PDGCodes.h:191
const int kPdgNuMu
Definition PDGCodes.h:30
@ kHadroSysGenErr
Definition GHepFlags.h:32