GENIEGenerator
Loading...
Searching...
No Matches
PhotonRESWdecPythia6.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#ifdef __GENIE_PYTHIA6_ENABLED__
14#if ROOT_VERSION_CODE >= ROOT_VERSION(5,15,6)
15#include <TMCParticle.h>
16#else
17#include <TMCParticle6.h>
18#endif
19#endif // __GENIE_PYTHIA6_ENABLED__
20
21//___________________________________________________________________________
23EventRecordVisitorI("genie::PhotonRESWdecPythia6")
24{
25 this->Initialize();
26 born = new Born();
27}
28//___________________________________________________________________________
30EventRecordVisitorI("genie::PhotonRESWdecPythia6", config)
31{
32 this->Initialize();
33 born = new Born();
34}
35//___________________________________________________________________________
40//____________________________________________________________________________
42{
43
44 if(!this->Wdecay(event)) {
45 LOG("PhotonRESGenerator", pWARN) << "W decayer failed!";
46 event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
48 exception.SetReason("Could not simulate the W decay system");
49 exception.SwitchOnFastForward();
50 throw exception;
51 return;
52 }
53
54}
55//___________________________________________________________________________
57#ifdef __GENIE_PYTHIA6_ENABLED__
58 event // avoid unused variable warning if PYTHIA6 is not enabled
59#endif
60) const
61{
62
63#ifdef __GENIE_PYTHIA6_ENABLED__
64 Interaction * interaction = event->Summary();
65 const InitialState & init_state = interaction->InitState();
66
67 //incoming v & struck particle & remnant nucleus
68 GHepParticle * nu = event->Probe();
69 GHepParticle * el = event->HitNucleon();
70
71 GHepParticle * target = event -> TargetNucleus();
72 if(target) event->AddParticle(target->Pdg(), kIStFinalStateNuclearRemnant, 1,-1,-1,-1, *(target->P4()), *(target->X4()) );
73
74 TVector3 unit_nu = nu->P4()->Vect().Unit();
75
76 int probepdg = init_state.ProbePdg();
77
78 long double Mtarget = init_state.Tgt().HitNucMass();
79 long double mlout = interaction->FSPrimLepton()->Mass(); //mass of charged lepton
80
81 long double Enuin = init_state.ProbeE(kRfLab);
82 long double s = born->GetS(Mtarget,Enuin);
83
84 long double n1 = interaction->Kine().GetKV(kKVn1);
85 long double n2 = interaction->Kine().GetKV(kKVn2);
86
87 long double costhCM = n1;
88 long double sinthCM = sqrtl(1-costhCM*costhCM);
89
90 long double xmin = fQ2PDFmin/2./Enuin/Mtarget;
91 long double x = expl( logl(xmin) + (logl(1.0)-logl(xmin))*n2 );
92 long double s_r = s*x;
93
94 // Boost velocity CM -> LAB
95 long double EnuinCM = sqrtl(s_r)/2.;
96 long double beta = (powl(Enuin,2)-powl(EnuinCM,2))/(powl(Enuin,2)+powl(EnuinCM,2));
97
98 // Final state primary lepton PDG code
99 int pdgl = interaction->FSPrimLeptonPdg();
100 assert(pdgl!=0);
101
102 if ( pdg::IsElectron(TMath::Abs(pdgl)) || pdg::IsMuon(TMath::Abs(pdgl)) || pdg::IsTau(TMath::Abs(pdgl)) ) {
103
104 long double ElpoutCM = (s_r+mlout*mlout)/sqrtl(s_r)/2.;
105 long double EnuoutCM = (s_r-mlout*mlout)/sqrtl(s_r)/2.;
106 LongLorentzVector p4_lpout( 0., EnuoutCM*sinthCM, EnuoutCM*costhCM, ElpoutCM );
107 LongLorentzVector p4_nuout( 0., -EnuoutCM*sinthCM, -EnuoutCM*costhCM, EnuoutCM );
108
109 p4_lpout.BoostZ(beta);
110 p4_nuout.BoostZ(beta);
111
112 TLorentzVector p4lp_o( (double)p4_lpout.Px(), (double)p4_lpout.Py(), (double)p4_lpout.Pz(), (double)p4_lpout.E() );
113 TLorentzVector p4nu_o( (double)p4_nuout.Px(), (double)p4_nuout.Py(), (double)p4_nuout.Pz(), (double)p4_nuout.E() );
114
115 // Randomize transverse components
117 double phi = 2* kPi * rnd->RndLep().Rndm();
118 p4lp_o.RotateZ(phi);
119 p4nu_o.RotateZ(phi);
120
121 //rotate from LAB=[0,0,Ev,Ev]->[px,py,pz,E]
122 p4lp_o.RotateUz(unit_nu);
123 p4nu_o.RotateUz(unit_nu);
124
125 int pdgvout = 0;
126 if ( pdg::IsElectron(pdgl) ) pdgvout = kPdgAntiNuE;
127 else if ( pdg::IsPositron(pdgl) ) pdgvout = kPdgNuE;
128 else if ( pdg::IsMuon(pdgl) ) pdgvout = kPdgAntiNuMu;
129 else if ( pdg::IsAntiMuon(pdgl) ) pdgvout = kPdgNuMu;
130 else if ( pdg::IsTau(pdgl) ) pdgvout = kPdgAntiNuTau;
131 else if ( pdg::IsAntiTau(pdgl) ) pdgvout = kPdgNuTau;
132
133 int pdgboson = pdg::IsNeutrino(probepdg) ? kPdgWP : kPdgWM;
134
135 // Create a GHepParticle and add it to the event record
136 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]
137 event->AddParticle( pdgl, kIStStableFinalState, 4, -1, -1, -1, p4lp_o, *(nu->X4()) );
138 event->AddParticle( pdgvout, kIStStableFinalState, 4, -1, -1, -1, p4nu_o, *(nu->X4()) );
139 event->Summary()->KinePtr()->SetFSLeptonP4(p4lp_o);
140
141 }
142 else {
143
144 char p6frame[10];
145 strcpy(p6frame, "CMS" );
146
147 char p6nu[10], p6tgt[10];
148 if (pdg::IsNeutrino(probepdg)) { strcpy(p6nu, "nu_e"); strcpy(p6tgt, "e+"); }
149 else { strcpy(p6nu, "nu_ebar"); strcpy(p6tgt, "e-"); }
150
151 int def61 = fPythia->GetMSTP(61);
152 int def71 = fPythia->GetMSTP(71);
153 int def206 = fPythia->GetMDME(206,1);
154 int def207 = fPythia->GetMDME(207,1);
155 int def208 = fPythia->GetMDME(208,1);
156 fPythia->SetMSTP(61,0); // (Default=2) master switch for initial-state QCD and QED radiation.
157 fPythia->SetMSTP(71,0); // (Default=2) master switch for initial-state QCD and QED radiation.
158 fPythia->SetMDME(206,1,0); //swicht off W decay leptonic modes
159 fPythia->SetMDME(207,1,0);
160 fPythia->SetMDME(208,1,0);
161
162 fPythia->Pyinit(p6frame, p6nu, p6tgt, sqrtl(s_r));
163 fPythia->Pyevnt();
164
165 fPythia->SetMSTP(61,def61);
166 fPythia->SetMSTP(71,def71);
167 fPythia->SetMDME(206,1,def206);
168 fPythia->SetMDME(207,1,def207);
169 fPythia->SetMDME(208,1,def208);
170
171 // get LUJETS record
172 fPythia->GetPrimaries();
173 TClonesArray * pythia_particles = (TClonesArray *) fPythia->ImportParticles("All");
174 int np = pythia_particles->GetEntries();
175 assert(np>0);
176
177 TMCParticle * particle = 0;
178 TIter piter(pythia_particles);
179 while( (particle = (TMCParticle *) piter.Next()) ) {
180
181 int pdgc = particle->GetKF();
182 int ks = particle->GetKS();
183
184 if ( ks==21 ) { continue; } //we dont want to save first particles from pythia (init states)
185
186 LongLorentzVector p4longo(particle->GetPx(), particle->GetPy(), particle->GetPz(), particle->GetEnergy());
187 p4longo.BoostZ(beta);
188
189 TLorentzVector p4o( (double)p4longo.Px(), (double)p4longo.Py(), (double)p4longo.Pz(), (double)p4longo.E() );
190 p4o.RotateUz(unit_nu);
191
192 TParticlePDG * part = PDGLibrary::Instance()->Find(pdgc);
193 if ( (ks==1 || ks==4) && p4o.E() < part->Mass() ) {
194 LOG("PhotonRESGenerator", pWARN) << "Putting at rest one stable particle generated by PYTHIA because E < m";
195 LOG("PhotonRESGenerator", pWARN) << "PDG = " << pdgc << " // State = " << ks;
196 LOG("PhotonRESGenerator", pWARN) << "E = " << p4o.E() << " // |p| = " << TMath::Sqrt(p4o.P());
197 LOG("PhotonRESGenerator", pWARN) << "p = [ " << p4o.Px() << " , " << p4o.Py() << " , " << p4o.Pz() << " ]";
198 LOG("PhotonRESGenerator", pWARN) << "m = " << p4o.M() << " // mpdg = " << part->Mass();
199 p4o.SetXYZT(0,0,0,part->Mass());
200 }
201
202 // copy final state particles to the event record
204
205 // fix numbering scheme used for mother/daughter assignments
206 int firstmother = -1;
207 int lastmother = -1;
208 int firstchild = -1;
209 int lastchild = -1;
210
211 if ( particle->GetParent() < 10 ) {
212 if ( TMath::Abs(pdgc)<7 ) { //outgoing quarks: mother will be the boson (saved in position 4)
213 firstmother = 4;
214 firstchild = particle->GetFirstChild() - 6;
215 lastchild = particle->GetLastChild() - 6;
216 }
217 else if ( TMath::Abs(pdgc)==24 ) { //produced W boson: mother will be the incoming neutrino
218 firstmother = 0;
219 firstchild = particle->GetFirstChild() - 6;
220 lastchild = particle->GetLastChild() - 6;
221 }
222 else if ( pdgc==22 ) { //radiative photons: mother will be the incoming electron
223 firstmother = 2;
224 }
225 }
226 else { //rest
227 firstmother = particle->GetParent() - 6; //shift to match boson position
228 firstchild = (particle->GetFirstChild()==0) ? particle->GetFirstChild() - 1 : particle->GetFirstChild() - 6;
229 lastchild = (particle->GetLastChild()==0) ? particle->GetLastChild() - 1 : particle->GetLastChild() - 6;
230 }
231
232 double lightspeed = 299792458e3; //c in mm/s. Used for time in PYTHIA t[s]=t_pythia[mm]/c[mm/s]
233 double vx = nu->X4()->X() + particle->GetVx()*1e12; //pythia gives position in [mm] while genie uses [fm]
234 double vy = nu->X4()->Y() + particle->GetVy()*1e12;
235 double vz = nu->X4()->Z() + particle->GetVz()*1e12;
236 double vt = nu->X4()->T() + particle->GetTime()/lightspeed;
237 TLorentzVector pos( vx, vy, vz, vt );
238
239 event->AddParticle(pdgc, ist, firstmother, lastmother, firstchild, lastchild, p4o, pos );
240
241 }
242
243 delete particle;
244 pythia_particles->Clear("C");
245
246 }
247
248 return true;
249#else
250 return false;
251#endif
252
253
254}
255//___________________________________________________________________________
261//____________________________________________________________________________
267//____________________________________________________________________________
269{
270
271#ifdef __GENIE_PYTHIA6_ENABLED__
272 GetParam( "SSBarSuppression", fSSBarSuppression );
273 GetParam( "GaussianPt2", fGaussianPt2 );
274 GetParam( "NonGaussianPt2Tail", fNonGaussianPt2Tail );
275 GetParam( "RemainingEnergyCutoff", fRemainingECutoff );
276 GetParam( "DiQuarkSuppression", fDiQuarkSuppression );
277 GetParam( "LightVMesonSuppression", fLightVMesonSuppression );
278 GetParam( "SVMesonSuppression", fSVMesonSuppression );
279 GetParam( "Lunda", fLunda );
280 GetParam( "Lundb", fLundb );
281 GetParam( "LundaDiq", fLundaDiq );
282
283 GetParam( "Q2Grid-Min", fQ2PDFmin );
284
285 // PYTHIA parameters only valid for HEDIS
286 double wmin; GetParam( "Xsec-Wmin", wmin ) ;
287 int warnings; GetParam( "Warnings", warnings ) ;
288 int errors; GetParam( "Errors", errors ) ;
289 int qrk_mass; GetParam( "QuarkMass", qrk_mass ) ;
290 fPythia->SetPARP(2, wmin); //(D = 10. GeV) lowest c.m. energy for the event as a whole that the program will accept to simulate. (bellow 2GeV pythia crashes)
291 fPythia->SetMSTU(26, warnings); // (Default=10) maximum number of warnings that are printed
292 fPythia->SetMSTU(22, errors); // (Default=10) maximum number of errors that are printed
293 fPythia->SetMSTJ(93, qrk_mass); // light (d, u, s, c, b) quark masses are taken from PARF(101) - PARF(105) rather than PMAS(1,1) - PMAS(5,1). Diquark masses are given as sum of quark masses, without spin splitting term.
294 fPythia->SetPMAS(24,1,kMw); //mass of the W boson (pythia=80.450 // genie=80.385)
295 fPythia->SetPMAS(24,2,0.); //set to 0 the width of the W boson to avoid problems with energy conservation
296 fPythia->SetPMAS(6,2,0.); //set to 0 the width of the top to avoid problems with energy conservation
297 fPythia->SetMDME(192,1,0); // W->dbar+t decay off
298 fPythia->SetMDME(196,1,0); // W->cbar+t decay off
299 fPythia->SetMDME(200,1,0); // W->cbar+t decay off
300
301 // PYTHIA tuned parameters
302 fPythia->SetPARJ(2, fSSBarSuppression );
303 fPythia->SetPARJ(21, fGaussianPt2 );
304 fPythia->SetPARJ(23, fNonGaussianPt2Tail );
305 fPythia->SetPARJ(33, fRemainingECutoff );
306 fPythia->SetPARJ(1, fDiQuarkSuppression );
307 fPythia->SetPARJ(11, fLightVMesonSuppression );
308 fPythia->SetPARJ(12, fSVMesonSuppression );
309 fPythia->SetPARJ(41, fLunda );
310 fPythia->SetPARJ(42, fLundb );
311 fPythia->SetPARJ(45, fLundaDiq );
312#endif // __GENIE_PYTHIA6_ENABLED__
313
314}
315//____________________________________________________________________________
317{
318
319#ifdef __GENIE_PYTHIA6_ENABLED__
320 fPythia = TPythia6::Instance();
321
322 // sync GENIE/PYTHIA6 seed number
324#endif
325
326}
#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
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
Initial State information.
const Target & Tgt(void) const
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)
double fNonGaussianPt2Tail
non gaussian pt2 tail parameterization
double fLundb
Lund b parameter.
double fLunda
Lund a parameter.
double fLundaDiq
adjustment of Lund a for di-quark
bool Wdecay(GHepRecord *event) const
double fSVMesonSuppression
strange vector meson suppression
double fSSBarSuppression
ssbar suppression
void ProcessEventRecord(GHepRecord *event) const
double fGaussianPt2
gaussian pt2 distribution width
double fLightVMesonSuppression
light vector meson suppression
void Configure(const Registry &config)
double fDiQuarkSuppression
di-quark suppression parameter
double fRemainingECutoff
remaining E cutoff stopping fragmentation
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
double HitNucMass(void) const
Definition Target.cxx:233
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
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
Simple functions for loading and reading nucleus dependent keys from config files.
@ kIStFinalStateNuclearRemnant
Definition GHepStatus.h:38
@ kIStDISPreFragmHadronicState
Definition GHepStatus.h:35
@ kIStStableFinalState
Definition GHepStatus.h:30
@ kIStDecayedState
Definition GHepStatus.h:32
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