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

PhotonRES W decay with pythia8. More...

#include <PhotonRESWdecPythia8.h>

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

Public Member Functions

 PhotonRESWdecPythia8 ()
 PhotonRESWdecPythia8 (string config)
 ~PhotonRESWdecPythia8 ()
void ProcessEventRecord (GHepRecord *event) const
void Configure (const Registry &config)
void Configure (string config)
Public Member Functions inherited from genie::EventRecordVisitorI
virtual ~EventRecordVisitorI ()
Public Member Functions inherited from genie::Algorithm
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.

Private Member Functions

bool Wdecay (GHepRecord *event) const
void Initialize (void) const
void LoadConfig (void)

Private Attributes

double fSSBarSuppression
 ssbar suppression
double fGaussianPt2
 gaussian pt2 distribution width
double fNonGaussianPt2Tail
 non gaussian pt2 tail parameterization
double fRemainingECutoff
 remaining E cutoff stopping fragmentation
double fDiQuarkSuppression
 di-quark suppression parameter
double fLightVMesonSuppression
 light vector meson suppression
double fSVMesonSuppression
 strange vector meson suppression
double fLunda
 Lund a parameter.
double fLundb
 Lund b parameter.
double fLundaDiq
 adjustment of Lund a for di-quark
double fQ2PDFmin
Bornborn

Additional Inherited Members

Static Public 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 Member Functions inherited from genie::EventRecordVisitorI
 EventRecordVisitorI ()
 EventRecordVisitorI (string name)
 EventRecordVisitorI (string name, string config)
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.
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

PhotonRES W decay with pythia8.

Author
Alfonso Garcia <aagarciasoto \at km3net.de> IFIC (Valencia)
Created:\n Dec 12, 2024
License:\n Copyright (c) 2003-2025, The GENIE Collaboration
For the full text of the license visit http://copyright.genie-mc.org or see $GENIE/LICENSE

Definition at line 46 of file PhotonRESWdecPythia8.h.

Constructor & Destructor Documentation

◆ PhotonRESWdecPythia8() [1/2]

PhotonRESWdecPythia8::PhotonRESWdecPythia8 ( )

Definition at line 14 of file PhotonRESWdecPythia8.cxx.

14 :
15EventRecordVisitorI("genie::PhotonRESWdecPythia8")
16{
17 this->Initialize();
18 born = new Born();
19}

References born, genie::EventRecordVisitorI::EventRecordVisitorI(), and Initialize().

◆ PhotonRESWdecPythia8() [2/2]

PhotonRESWdecPythia8::PhotonRESWdecPythia8 ( string config)

Definition at line 21 of file PhotonRESWdecPythia8.cxx.

21 :
22EventRecordVisitorI("genie::PhotonRESWdecPythia8", config)
23{
24 this->Initialize();
25 born = new Born();
26}

References born, genie::EventRecordVisitorI::EventRecordVisitorI(), and Initialize().

◆ ~PhotonRESWdecPythia8()

PhotonRESWdecPythia8::~PhotonRESWdecPythia8 ( )

Definition at line 28 of file PhotonRESWdecPythia8.cxx.

29{
30
31}

Member Function Documentation

◆ Configure() [1/2]

void PhotonRESWdecPythia8::Configure ( const Registry & config)
virtual

Configure the algorithm with an external registry The registry is merged with the top level registry if it is owned, Otherwise a copy of it is added with the highest priority

Reimplemented from genie::Algorithm.

Definition at line 225 of file PhotonRESWdecPythia8.cxx.

226{
227 Algorithm::Configure(config);
228 this->LoadConfig();
229}
virtual void Configure(const Registry &config)
Definition Algorithm.cxx:62

References genie::Algorithm::Configure(), and LoadConfig().

◆ Configure() [2/2]

void PhotonRESWdecPythia8::Configure ( string config)
virtual

Configure the algorithm from the AlgoConfigPool based on param_set string given in input An algorithm contains a vector of registries coming from different xml configuration files, which are loaded according a very precise prioriy This methods will load a number registries in order of priority: 1) "Tunable" parameter set from CommonParametes. This is loaded with the highest prioriry and it is designed to be used for tuning procedure Usage not expected from the user. 2) For every string defined in "CommonParame" the corresponding parameter set will be loaded from CommonParameter.xml 3) parameter set specified by the config string and defined in the xml file of the algorithm 4) if config is not "Default" also the Default parameter set from the same xml file will be loaded Effectively this avoids the repetion of a parameter when it is not changed in the requested configuration

Reimplemented from genie::Algorithm.

Definition at line 231 of file PhotonRESWdecPythia8.cxx.

232{
233 Algorithm::Configure(config);
234 this->LoadConfig();
235}

References genie::Algorithm::Configure(), and LoadConfig().

◆ Initialize()

void PhotonRESWdecPythia8::Initialize ( void ) const
private

Definition at line 350 of file PhotonRESWdecPythia8.cxx.

351{
352
353#ifdef __GENIE_PYTHIA8_ENABLED__
354 //One cool feature of having independent Pythia8 instances
355 //We can define our intial state with the proper setting (i.e. disabling decays)
356 //without affecting other classes of GENIE that also use Pythia8
357
358 fPythiaP = new Pythia8::Pythia();
359 fPythiaP->readString("Print:quiet = on");
360 fPythiaP->readString("Random:setSeed = on");
361 fPythiaP->readString("WeakSingleBoson:ffbar2ffbar(s:W) = on");
362 fPythiaP->readString("PDF:lepton = off");
363 fPythiaP->readString("24:onMode = off");
364 fPythiaP->readString("24:onIfAny = 1 2 3 4 5"); //enable W->hadron only
365 fPythiaP->readString("Beams:idA = 12");
366 fPythiaP->readString("Beams:idB = -11");
367
368 fPythiaN = new Pythia8::Pythia();
369 fPythiaN->readString("Print:quiet = on");
370 fPythiaN->readString("Random:setSeed = on");
371 fPythiaN->readString("WeakSingleBoson:ffbar2ffbar(s:W) = on");
372 fPythiaN->readString("PDF:lepton = off");
373 fPythiaN->readString("24:onMode = off");
374 fPythiaN->readString("24:onIfAny = 1 2 3 4 5"); //enable W->hadron only
375 fPythiaN->readString("Beams:idA = -12");
376 fPythiaN->readString("Beams:idB = 11");
377#endif
378
379}

Referenced by PhotonRESWdecPythia8(), and PhotonRESWdecPythia8().

◆ LoadConfig()

void PhotonRESWdecPythia8::LoadConfig ( void )
private

Definition at line 237 of file PhotonRESWdecPythia8.cxx.

238{
239
240#ifdef __GENIE_PYTHIA8_ENABLED__
241 GetParam( "SSBarSuppression", fSSBarSuppression );
242 GetParam( "GaussianPt2", fGaussianPt2 );
243 GetParam( "NonGaussianPt2Tail", fNonGaussianPt2Tail );
244 GetParam( "RemainingEnergyCutoff", fRemainingECutoff );
245 GetParam( "DiQuarkSuppression", fDiQuarkSuppression );
246 GetParam( "LightVMesonSuppression", fLightVMesonSuppression );
247 GetParam( "SVMesonSuppression", fSVMesonSuppression );
248 GetParam( "Lunda", fLunda );
249 GetParam( "Lundb", fLundb );
250 GetParam( "LundaDiq", fLundaDiq );
251
252 GetParam( "Q2Grid-Min", fQ2PDFmin );
253
254 fPythiaP->settings.parm("StringFlav:probStoUD", fSSBarSuppression);
255 fPythiaP->settings.parm("Diffraction:primKTwidth", fGaussianPt2);
256 fPythiaP->settings.parm("StringPT:enhancedFraction", fNonGaussianPt2Tail);
257 fPythiaP->settings.parm("StringFragmentation:stopMass", fRemainingECutoff);
258 fPythiaP->settings.parm("StringFlav:probQQtoQ", fDiQuarkSuppression);
259 fPythiaP->settings.parm("StringFlav:mesonUDvector", fLightVMesonSuppression);
260 fPythiaP->settings.parm("StringFlav:mesonSvector", fSVMesonSuppression);
261 fPythiaP->settings.parm("StringZ:aLund", fLunda);
262 fPythiaP->settings.parm("StringZ:bLund", fLundb);
263 fPythiaP->settings.parm("StringZ:aExtraDiquark", fLundaDiq);
264
265 fPythiaN->settings.parm("StringFlav:probStoUD", fSSBarSuppression);
266 fPythiaN->settings.parm("Diffraction:primKTwidth", fGaussianPt2);
267 fPythiaN->settings.parm("StringPT:enhancedFraction", fNonGaussianPt2Tail);
268 fPythiaN->settings.parm("StringFragmentation:stopMass", fRemainingECutoff);
269 fPythiaN->settings.parm("StringFlav:probQQtoQ", fDiQuarkSuppression);
270 fPythiaN->settings.parm("StringFlav:mesonUDvector", fLightVMesonSuppression);
271 fPythiaN->settings.parm("StringFlav:mesonSvector", fSVMesonSuppression);
272 fPythiaN->settings.parm("StringZ:aLund", fLunda);
273 fPythiaN->settings.parm("StringZ:bLund", fLundb);
274 fPythiaN->settings.parm("StringZ:aExtraDiquark", fLundaDiq);
275
276 // Same default mass of the W boson in pythia8 and genie, so no need to change in pythia8
277 // No problem with energy conservation W and top decays, so no need to set the width to 0
278
279 // Important for not decaying too early!!!
280 // Copy of the method from Decayer class to have a consistent treatment.
281 // The idea is to not decay in this stage what we dont want to decay in Decayer.
282 // Hadronization -> Decayer
283 //------------------------------
284 // Decay -> Decay
285 // Undecay -> Undecay
286 // Undecay -> Decay
287 // Decay -> Undecay -> This is the problematic one
288 // So we loop over all the particles for which we inhibit decay in the config file
289 // for the Decayer stage and we inhibit it also at hadronization.
290 // This is not need in LeptonHadronization and PhotonCOH because they use the
291 // Pythia8Singleton class which already defines the decays consistently
292 //
293 RgKeyList klist = GetConfig().FindKeys("DecayParticleWithCode=");
294 RgKeyList::const_iterator kiter = klist.begin();
295 for( ; kiter != klist.end(); ++kiter) {
296 RgKey key = *kiter;
297 bool decay = GetConfig().GetBool(key);
298 vector<string> kv = utils::str::Split(key,"=");
299 int pdg_code = atoi(kv[1].c_str());
300 if(!decay) {
301 LOG("PhotonRESGenerator", pDEBUG) << "Configured to inhibit decays for " << pdg_code;
302 auto pdentryP = fPythiaP->particleData.particleDataEntryPtr(pdg_code);
303 for (int ichan=0; ichan<=pdentryP->sizeChannels()-1; ++ichan) {
304 int onMode = pdentryP->channel(ichan).onMode();
305 bool is_particle = (pdg_code>0);
306 switch ( onMode ) {
307 case 0: // off for particles & antiparticles
308 break; // already off
309 case 1: // on for both particle & antiparticle
310 onMode = (is_particle ? 3 : 2);
311 break;
312 case 2: // on for particle only
313 onMode = (is_particle ? 0 : 2);
314 break;
315 case 3: // on for antiparticle only
316 onMode = (is_particle ? 3 : 0);
317 break;
318 }
319 pdentryP->channel(ichan).onMode(onMode);
320 }
321
322 auto pdentryN = fPythiaN->particleData.particleDataEntryPtr(pdg_code);
323 for (int ichan=0; ichan<=pdentryN->sizeChannels()-1; ++ichan) {
324 int onMode = pdentryN->channel(ichan).onMode();
325 bool is_particle = (pdg_code>0);
326 switch ( onMode ) {
327 case 0: // off for particles & antiparticles
328 break; // already off
329 case 1: // on for both particle & antiparticle
330 onMode = (is_particle ? 3 : 2);
331 break;
332 case 2: // on for particle only
333 onMode = (is_particle ? 0 : 2);
334 break;
335 case 3: // on for antiparticle only
336 onMode = (is_particle ? 3 : 0);
337 break;
338 }
339 pdentryN->channel(ichan).onMode(onMode);
340 }
341
342 }
343 }
344
345 //Important: We dont initialize Pythia8 here because it must be done event by event. See above
346#endif
347
348}
#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
string RgKey
virtual const Registry & GetConfig(void) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
double fLundaDiq
adjustment of Lund a for di-quark
double fNonGaussianPt2Tail
non gaussian pt2 tail parameterization
double fSSBarSuppression
ssbar suppression
double fGaussianPt2
gaussian pt2 distribution width
double fLundb
Lund b parameter.
double fLightVMesonSuppression
light vector meson suppression
double fSVMesonSuppression
strange vector meson suppression
double fDiQuarkSuppression
di-quark suppression parameter
double fRemainingECutoff
remaining E cutoff stopping fragmentation
double fLunda
Lund a parameter.
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
vector< string > Split(string input, string delim)
vector< RgKey > RgKeyList
Definition Registry.h:50

References fDiQuarkSuppression, fGaussianPt2, genie::Registry::FindKeys(), fLightVMesonSuppression, fLunda, fLundaDiq, fLundb, fNonGaussianPt2Tail, fQ2PDFmin, fRemainingECutoff, fSSBarSuppression, fSVMesonSuppression, genie::Registry::GetBool(), genie::Algorithm::GetConfig(), genie::Algorithm::GetParam(), LOG, pDEBUG, and genie::utils::str::Split().

Referenced by Configure(), and Configure().

◆ ProcessEventRecord()

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

Implements genie::EventRecordVisitorI.

Definition at line 33 of file PhotonRESWdecPythia8.cxx.

34{
35
36 if(!this->Wdecay(event)) {
37 LOG("PhotonRESGenerator", pWARN) << "W decayer failed!";
38 event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
39 genie::exceptions::EVGThreadException exception;
40 exception.SetReason("Could not simulate the W decay system");
41 exception.SwitchOnFastForward();
42 throw exception;
43 return;
44 }
45
46}
#define pWARN
Definition Messenger.h:60
bool Wdecay(GHepRecord *event) const
@ kHadroSysGenErr
Definition GHepFlags.h:32

References genie::kHadroSysGenErr, LOG, pWARN, genie::exceptions::EVGThreadException::SetReason(), genie::exceptions::EVGThreadException::SwitchOnFastForward(), and Wdecay().

◆ Wdecay()

bool PhotonRESWdecPythia8::Wdecay ( GHepRecord * event) const
private

Definition at line 48 of file PhotonRESWdecPythia8.cxx.

53{
54
55#ifdef __GENIE_PYTHIA8_ENABLED__
56 Interaction * interaction = event->Summary();
57 const InitialState & init_state = interaction->InitState();
58
59 //incoming v & struck particle & remnant nucleus
60 GHepParticle * nu = event->Probe();
61 GHepParticle * el = event->HitNucleon();
62
63 GHepParticle * target = event -> TargetNucleus();
64 if(target) event->AddParticle(target->Pdg(), kIStFinalStateNuclearRemnant, 1,-1,-1,-1, *(target->P4()), *(target->X4()) );
65
66 TVector3 unit_nu = nu->P4()->Vect().Unit();
67
68 int probepdg = init_state.ProbePdg();
69
70 long double Mtarget = init_state.Tgt().HitNucMass();
71 long double mlout = interaction->FSPrimLepton()->Mass(); //mass of charged lepton
72
73 long double Enuin = init_state.ProbeE(kRfLab);
74 long double s = born->GetS(Mtarget,Enuin);
75
76 long double n1 = interaction->Kine().GetKV(kKVn1);
77 long double n2 = interaction->Kine().GetKV(kKVn2);
78
79 long double costhCM = n1;
80 long double sinthCM = sqrtl(1-costhCM*costhCM);
81
82 long double xmin = fQ2PDFmin/2./Enuin/Mtarget;
83 long double x = expl( logl(xmin) + (logl(1.0)-logl(xmin))*n2 );
84 long double s_r = s*x;
85
86 // Boost velocity CM -> LAB
87 long double EnuinCM = 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
108 RandomGen * rnd = RandomGen::Instance();
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(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
137 //We have to initialize Pythia8 everytime we generate a new event
138 //And we have to change the seed because it would regenerate the same event
139 //if we dont do it.
140 RandomGen * rnd = RandomGen::Instance();
141 Pythia8::Event fEvent;
142 if (pdg::IsNeutrino(probepdg)) {
143 fPythiaP->settings.mode("Random:seed", rnd->RndLep().Integer(100000000));
144 fPythiaP->settings.parm("Beams:eCM",sqrtl(s_r));
145 fPythiaP->init();
146 fPythiaP->next();
147 // fPythiaP->event.list();
148 // fPythiaP->stat();
149 fEvent = fPythiaP->event;
150 }
151 else {
152 fPythiaN->settings.mode("Random:seed", rnd->RndLep().Integer(100000000));
153 fPythiaN->settings.parm("Beams:eCM",sqrtl(s_r));
154 fPythiaN->init();
155 fPythiaN->next();
156 // fPythiaN->event.list();
157 // fPythiaN->stat();
158 fEvent = fPythiaN->event;
159 }
160
161 int np = fEvent.size();
162 assert(np>0);
163
164 for (int i = 5; i < np; i++) { // ignore first entry -> system + input particles
165
166 int pdgc = fEvent[i].id();
167 if (!PDGLibrary::Instance()->Find(pdgc)) continue; // some intermediate particles not part of genie tables
168
169 int ks = fEvent[i].status();
170
171 LongLorentzVector p4longo(fEvent[i].px(), fEvent[i].py(), fEvent[i].pz(), fEvent[i].e());
172 p4longo.BoostZ(beta);
173
174 TLorentzVector p4o( (double)p4longo.Px(), (double)p4longo.Py(), (double)p4longo.Pz(), (double)p4longo.E() );
175 p4o.RotateUz(unit_nu);
176
177 double massPDG = PDGLibrary::Instance()->Find(pdgc)->Mass();
178 if ( ks>0 && p4o.E()<massPDG ) {
179 LOG("PhotonRESGenerator", pWARN) << "Putting at rest one stable particle generated by PYTHIA because E < m";
180 LOG("PhotonRESGenerator", pWARN) << "PDG = " << pdgc << " // State = " << ks;
181 LOG("PhotonRESGenerator", pWARN) << "E = " << p4o.E() << " // |p| = " << TMath::Sqrt(p4o.P());
182 LOG("PhotonRESGenerator", pWARN) << "p = [ " << p4o.Px() << " , " << p4o.Py() << " , " << p4o.Pz() << " ]";
183 LOG("PhotonRESGenerator", pWARN) << "m = " << p4o.M() << " // mpdg = " << massPDG;
184 p4o.SetXYZT(0,0,0,massPDG);
185 }
186
187 // copy final state particles to the event record
189
190 // fix numbering scheme used for mother/daughter assignments
191 int firstmother = -1;
192 int lastmother = -1;
193 int firstchild = -1;
194 int lastchild = -1;
195
196 if ( fEvent[i].mother1() == 3 ) { //produced W boson: mother will be the incoming neutrino
197 firstmother = 0;
198 firstchild = fEvent[i].daughter1() - 3;
199 lastchild = fEvent[i].daughter2() - 3;
200 }
201 else { //rest
202 firstmother = fEvent[i].mother1() - 3; //shift to match boson position
203 firstchild = fEvent[i].daughter1() - 3;
204 lastchild = fEvent[i].daughter2() - 3;
205 }
206
207 double vx = nu->X4()->X() + fEvent[i].xProd()*1e12; //pythia gives position in [mm] while genie uses [fm]
208 double vy = nu->X4()->Y() + fEvent[i].yProd()*1e12;
209 double vz = nu->X4()->Z() + fEvent[i].zProd()*1e12;
210 double vt = nu->X4()->T() + fEvent[i].tProd()*(units::millimeter/units::second);
211 TLorentzVector pos( vx, vy, vz, vt );
212
213 event->AddParticle(pdgc, ist, firstmother, lastmother, firstchild, lastchild, p4o, pos );
214 }
215
216 }
217
218 return true;
219#else
220 return false;
221#endif
222
223}
int Pdg(void) const
const TLorentzVector * P4(void) const
const TLorentzVector * X4(void) const
const Target & Tgt(void) const
TParticlePDG * Probe(void) const
int ProbePdg(void) const
double ProbeE(RefFrame_t rf) const
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)
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
double HitNucMass(void) const
Definition Target.cxx:233
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
@ 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

References genie::utils::math::LongLorentzVector::BoostZ(), born, genie::utils::math::LongLorentzVector::E(), e, genie::PDGLibrary::Find(), fQ2PDFmin, genie::Interaction::FSPrimLepton(), genie::Interaction::FSPrimLeptonPdg(), genie::Kinematics::GetKV(), genie::Target::HitNucMass(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::RandomGen::Instance(), genie::pdg::IsAntiMuon(), genie::pdg::IsAntiTau(), genie::pdg::IsElectron(), genie::pdg::IsMuon(), genie::pdg::IsNeutrino(), genie::pdg::IsPositron(), genie::pdg::IsTau(), genie::Interaction::Kine(), genie::kIStDecayedState, genie::kIStDISPreFragmHadronicState, genie::kIStFinalStateNuclearRemnant, genie::kIStStableFinalState, genie::kKVn1, genie::kKVn2, genie::kPdgAntiNuE, genie::kPdgAntiNuMu, genie::kPdgAntiNuTau, genie::kPdgNuE, genie::kPdgNuMu, genie::kPdgNuTau, genie::kPdgWM, genie::kPdgWP, genie::constants::kPi, genie::kRfLab, LOG, genie::units::millimeter, genie::GHepParticle::P4(), genie::GHepParticle::Pdg(), genie::InitialState::ProbeE(), genie::InitialState::ProbePdg(), pWARN, genie::utils::math::LongLorentzVector::Px(), genie::utils::math::LongLorentzVector::Py(), genie::utils::math::LongLorentzVector::Pz(), genie::RandomGen::RndLep(), genie::units::second, genie::InitialState::Tgt(), and genie::GHepParticle::X4().

Referenced by ProcessEventRecord().

Member Data Documentation

◆ born

Born* genie::PhotonRESWdecPythia8::born
private

Definition at line 82 of file PhotonRESWdecPythia8.h.

Referenced by PhotonRESWdecPythia8(), PhotonRESWdecPythia8(), and Wdecay().

◆ fDiQuarkSuppression

double genie::PhotonRESWdecPythia8::fDiQuarkSuppression
private

di-quark suppression parameter

Definition at line 73 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fGaussianPt2

double genie::PhotonRESWdecPythia8::fGaussianPt2
private

gaussian pt2 distribution width

Definition at line 70 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fLightVMesonSuppression

double genie::PhotonRESWdecPythia8::fLightVMesonSuppression
private

light vector meson suppression

Definition at line 74 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fLunda

double genie::PhotonRESWdecPythia8::fLunda
private

Lund a parameter.

Definition at line 76 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fLundaDiq

double genie::PhotonRESWdecPythia8::fLundaDiq
private

adjustment of Lund a for di-quark

Definition at line 78 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fLundb

double genie::PhotonRESWdecPythia8::fLundb
private

Lund b parameter.

Definition at line 77 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fNonGaussianPt2Tail

double genie::PhotonRESWdecPythia8::fNonGaussianPt2Tail
private

non gaussian pt2 tail parameterization

Definition at line 71 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fQ2PDFmin

double genie::PhotonRESWdecPythia8::fQ2PDFmin
private

Definition at line 80 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig(), and Wdecay().

◆ fRemainingECutoff

double genie::PhotonRESWdecPythia8::fRemainingECutoff
private

remaining E cutoff stopping fragmentation

Definition at line 72 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fSSBarSuppression

double genie::PhotonRESWdecPythia8::fSSBarSuppression
private

ssbar suppression

Definition at line 69 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().

◆ fSVMesonSuppression

double genie::PhotonRESWdecPythia8::fSVMesonSuppression
private

strange vector meson suppression

Definition at line 75 of file PhotonRESWdecPythia8.h.

Referenced by LoadConfig().


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