Public Member Functions | Private Member Functions | Private Attributes | List of all members
genie::GLRESGenerator Class Reference

Glashow resonance event generator. More...

#include <GLRESGenerator.h>

Inheritance diagram for genie::GLRESGenerator:
genie::EventRecordVisitorI genie::Algorithm

Public Member Functions

 GLRESGenerator ()
 
 GLRESGenerator (string config)
 
 ~GLRESGenerator ()
 
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. More...
 
virtual AlgStatus_t GetStatus (void) const
 Get algorithm status. More...
 
virtual bool AllowReconfig (void) const
 
virtual AlgCmp_t Compare (const Algorithm *alg) const
 Compare with input algorithm. More...
 
virtual void SetId (const AlgId &id)
 Set algorithm ID. More...
 
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. More...
 

Private Member Functions

void LoadConfig (void)
 

Private Attributes

TPythia6 * fPythia
 PYTHIA6 wrapper class. More...
 
double fWmin
 

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)
 
- 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. More...
 
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. More...
 
int GetParamVectKeys (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 More...
 
int AddLowRegistry (Registry *rp, bool owns=true)
 add registry with lowest priority, also update ownership More...
 
int MergeTopRegistry (const Registry &r)
 
int AddTopRegisties (const vector< Registry * > &rs, bool owns=false)
 Add registries with top priority, also udated Ownerships. More...
 
- Protected Attributes inherited from genie::Algorithm
bool fAllowReconfig
 
bool fOwnsSubstruc
 true if it owns its substructure (sub-algs,...) More...
 
AlgId fID
 algorithm name and configuration set More...
 
vector< Registry * > fConfVect
 
vector< boolfOwnerships
 ownership for every registry in fConfVect More...
 
AlgStatus_t fStatus
 algorithm execution status More...
 
AlgMapfOwnedSubAlgMp
 local pool for owned sub-algs (taken out of the factory pool) More...
 

Detailed Description

Glashow resonance event generator.

Author
Costas Andreopoulos <constantinos.andreopoulos cern.ch> University of Liverpool & STFC Rutherford Appleton Laboratory

Feb 15, 2008

Copyright (c) 2003-2020, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 31 of file GLRESGenerator.h.

Constructor & Destructor Documentation

GLRESGenerator::GLRESGenerator ( )

Definition at line 44 of file GLRESGenerator.cxx.

44  :
45 EventRecordVisitorI("genie::GLRESGenerator")
46 {
47 
48 }
GLRESGenerator::GLRESGenerator ( string  config)

Definition at line 50 of file GLRESGenerator.cxx.

50  :
51 EventRecordVisitorI("genie::GLRESGenerator", config)
52 {
53 
54 }
static Config * config
Definition: config.cpp:1054
GLRESGenerator::~GLRESGenerator ( )

Definition at line 56 of file GLRESGenerator.cxx.

57 {
58 
59 }

Member Function Documentation

void GLRESGenerator::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 290 of file GLRESGenerator.cxx.

291 {
292  Algorithm::Configure(config);
293  this->LoadConfig();
294 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void GLRESGenerator::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 296 of file GLRESGenerator.cxx.

297 {
299  this->LoadConfig();
300 }
static Config * config
Definition: config.cpp:1054
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void GLRESGenerator::LoadConfig ( void  )
private

Definition at line 302 of file GLRESGenerator.cxx.

303 {
304 #ifdef __GENIE_PYTHIA6_ENABLED__
305  fPythia = TPythia6::Instance();
306 
307  // sync GENIE/PYTHIA6 seed number
309 
310  // PYTHIA parameters only valid for HEDIS
311  GetParam( "Xsec-Wmin", fWmin ) ;
312  fPythia->SetPARP(2, fWmin); //(D = 10. GeV) lowest c.m. energy for the event as a whole that the program will accept to simulate. (bellow 2GeV pythia crashes)
313 
314  int warnings; GetParam( "Warnings", warnings ) ;
315  int errors; GetParam( "Errors", errors ) ;
316  int qrk_mass; GetParam( "QuarkMass", qrk_mass ) ;
317  fPythia->SetMSTU(26, warnings); // (Default=10) maximum number of warnings that are printed
318  fPythia->SetMSTU(22, errors); // (Default=10) maximum number of errors that are printed
319  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.
320 
321 #endif // __GENIE_PYTHIA6_ENABLED__
322 
323 }
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
TPythia6 * fPythia
PYTHIA6 wrapper class.
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
void GLRESGenerator::ProcessEventRecord ( GHepRecord event) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 61 of file GLRESGenerator.cxx.

66 {
67 
68 #ifdef __GENIE_PYTHIA6_ENABLED__
69 
70  Interaction * interaction = event->Summary();
71  const InitialState & init_state = interaction->InitState();
72 
73  //incoming v & struck particle & remnant nucleus
74  GHepParticle * nu = event -> Probe();
75  GHepParticle * el = event -> HitElectron();
76  GHepParticle * target = event -> TargetNucleus();
77  assert(nu);
78  assert(el);
79  assert(target);
80 
81  if(target) event->AddParticle(target->Pdg(), kIStStableFinalState, 1,-1,-1,-1, *(target->P4()), *(target->X4()) );
82 
83  LongLorentzVector p4_nu( * event->Probe()->P4() );
84  LongLorentzVector p4_el( * event->HitElectron()->P4() );
85  LongLorentzVector p4_Wlong( p4_nu.Px()+p4_el.Px(), p4_nu.Py()+p4_el.Py(), p4_nu.Pz()+p4_el.Pz(), p4_nu.E()+p4_el.E() );
86 
87  long double Wmass = p4_Wlong.M();
88  LOG("GLRESGenerator", pINFO) << "Wmass : " << Wmass;
89 
90  if(Wmass < fWmin) {
91  LOG("GLRESGenerator", pWARN) << "Low invariant mass, W = " << Wmass << " GeV!!";
92  event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
94  exception.SetReason("Could not simulate the hadronic system");
95  exception.SwitchOnFastForward();
96  throw exception;
97  return;
98  }
99 
100  // Final state primary lepton PDG code
101  int pdgl = interaction->FSPrimLeptonPdg();
102  assert(pdgl!=0);
103 
104  LOG("GLRESGenerator", pINFO) << "Channel : " << interaction->FSPrimLeptonPdg();
105 
106  if ( pdg::IsElectron(pdgl) || pdg::IsMuon(pdgl) || pdg::IsTau(pdgl) ) {
107 
108  // Get selected kinematics
109  long double y = interaction->Kine().y(true);
110  assert(y>0 && y<1);
111 
112  // Compute the neutrino and muon energy
113  long double Ev = init_state.ProbeE(kRfLab);
114  long double El = y*Ev;
115 
116  LOG("GLRESGenerator", pINFO) << "Ev = " << Ev << ", y = " << y << ", -> El = " << El;
117 
118  // Compute the momentum transfer and scattering angle
119  long double El2 = powl(El,2);
120  long double ml = interaction->FSPrimLepton()->Mass();
121  long double ml2 = powl(ml,2);
122  long double pl = sqrtl(El2-ml2);
123 
124  long double Q2 = 2*(Ev-El)*kElectronMass + kElectronMass*kElectronMass;
125  long double costh = (El-0.5*(Q2+ml2)/Ev)/pl;
126  long double sinth = sqrtl( 1-powl(costh,2.) );
127  // Randomize transverse components
128  RandomGen * rnd = RandomGen::Instance();
129  long double phi = 2* constants::kPi * rnd->RndLep().Rndm();
130 
131  LOG("GLRESGenerator", pINFO) << "Q2 = " << Q2 << ", cos(theta) = " << costh << ", phi = " << phi;
132 
133  long double plx = pl*sinth*cosl(phi);
134  long double ply = pl*sinth*sinl(phi);
135  long double plz = pl*costh;
136 
137  // Lepton 4-momentum in the LAB frame
138  LongLorentzVector p4lplong( plx, ply, plz, El );
139  p4lplong.Rotate(p4_nu);
140  LongLorentzVector p4nulong( p4_Wlong.Px()-p4lplong.Px(), p4_Wlong.Py()-p4lplong.Py(), p4_Wlong.Pz()-p4lplong.Pz(), p4_Wlong.E()-p4lplong.E() );
141 
142  TLorentzVector p4_W ( (double)p4_Wlong.Px(), (double)p4_Wlong.Py(), (double)p4_Wlong.Pz(), (double)p4_Wlong.E() );
143  TLorentzVector p4_lpout( (double)p4lplong.Px(), (double)p4lplong.Py(), (double)p4lplong.Pz(), (double)p4lplong.E() );
144  TLorentzVector p4_nuout( (double)p4nulong.Px(), (double)p4nulong.Py(), (double)p4nulong.Pz(), (double)p4nulong.E() );
145 
146  int pdgvout = 0;
147  if ( pdg::IsElectron(pdgl) ) pdgvout = kPdgAntiNuE;
148  else if ( pdg::IsMuon(pdgl) ) pdgvout = kPdgAntiNuMu;
149  else if ( pdg::IsTau(pdgl) ) pdgvout = kPdgAntiNuTau;
150 
151  // Create a GHepParticle and add it to the event record
152  event->AddParticle( kPdgWM, kIStDecayedState, 0, -1, 5, 6, p4_W, *(nu->X4()) ); //W [mothers: nuebar_in,e_in][daugthers: nulbar_out,lp_out]
153  event->AddParticle( pdgvout, kIStStableFinalState, 4, -1, -1, -1, p4_nuout, *(nu->X4()) );
154  event->AddParticle( pdgl, kIStStableFinalState, 4, -1, -1, -1, p4_lpout, *(nu->X4()) );
155  event->Summary()->KinePtr()->SetFSLeptonP4(p4_lpout);
156 
157  }
158  else {
159 
160  char p6frame[10], p6nu[10], p6tgt[10];
161  strcpy(p6frame, "CMS" );
162  strcpy(p6nu, "nu_ebar");
163  strcpy(p6tgt, "e-" );
164 
165 
166  int mstp61 = fPythia->GetMSTP(61);
167  int mstp71 = fPythia->GetMSTP(71);
168  int mdme206 = fPythia->GetMDME(206,1);
169  int mdme207 = fPythia->GetMDME(207,1);
170  int mdme208 = fPythia->GetMDME(208,1);
171  double pmas1_W = fPythia->GetPMAS(24,1);
172  double pmas2_W = fPythia->GetPMAS(24,2);
173  double pmas2_t = fPythia->GetPMAS(6,2);
174  fPythia->SetMSTP(61,0); // (Default=2) master switch for initial-state QCD and QED radiation.
175  fPythia->SetMSTP(71,0); // (Default=2) master switch for initial-state QCD and QED radiation.
176  fPythia->SetMDME(206,1,0); //swicht off W decay leptonic modes
177  fPythia->SetMDME(207,1,0);
178  fPythia->SetMDME(208,1,0);
179  fPythia->SetPMAS(24,1,kMw); //mass of the W boson (pythia=80.450 // genie=80.385)
180  fPythia->SetPMAS(24,2,0.); //set to 0 the width of the W boson to avoid problems with energy conservation
181  fPythia->SetPMAS(6,2,0.); //set to 0 the width of the top to avoid problems with energy conservation
182 
183  fPythia->Pyinit(p6frame, p6nu, p6tgt, (double)Wmass);
184 
185  fPythia->Pyevnt();
186 
187  fPythia->SetMSTP(61,mstp61);
188  fPythia->SetMSTP(71,mstp71);
189  fPythia->SetMDME(206,1,mdme206);
190  fPythia->SetMDME(207,1,mdme207);
191  fPythia->SetMDME(208,1,mdme208);
192  fPythia->SetPMAS(24,1,pmas1_W);
193  fPythia->SetPMAS(24,2,pmas2_W);
194  fPythia->SetPMAS(6,2,pmas2_t);
195 
196 
197  // Use for debugging purposes
198  //fPythia->Pylist(3);
199 
200  // get LUJETS record
201  fPythia->GetPrimaries();
202  TClonesArray * pythia_particles = (TClonesArray *) fPythia->ImportParticles("All");
203  // copy PYTHIA container to a new TClonesArray so as to transfer ownership
204  // of the container and of its elements to the calling method
205  int np = pythia_particles->GetEntries();
206  assert(np>0);
207  TClonesArray * particle_list = new TClonesArray("genie::GHepParticle", np);
208  particle_list->SetOwner(true);
209 
210  // Boost velocity LAB' -> Resonance rest frame
211  long double beta = p4_Wlong.P()/p4_Wlong.E();
212 
213  TMCParticle * p = 0;
214  TIter particle_iter(pythia_particles);
215  while( (p = (TMCParticle *) particle_iter.Next()) ) {
216 
217  int pdgc = p->GetKF();
218  int ks = p->GetKS();
219  int parent = p->GetParent();
220 
221  if ( ks==21 ) { continue; } //we dont want to save first particles from pythia (init states)
222 
223  LongLorentzVector p4long( p->GetPx(), p->GetPy(), p->GetPz(), p->GetEnergy() );
224  p4long.Boost(beta);
225  p4long.Rotate(p4_Wlong);
226 
227  TLorentzVector p4( (double)p4long.Px(), (double)p4long.Py(), (double)p4long.Pz(), (double)p4long.E() );
228 
229  double massPDG = PDGLibrary::Instance()->Find(pdgc)->Mass();
230  if ( (ks==1 || ks==4) && p4.E() < massPDG ) {
231  LOG("GLRESGenerator", pWARN) << "Putting at rest one stable particle generated by PYTHIA because E < m";
232  LOG("GLRESGenerator", pWARN) << "PDG = " << pdgc << " // State = " << ks;
233  LOG("GLRESGenerator", pWARN) << "E = " << p4.E() << " // |p| = " << TMath::Sqrt(p4.P());
234  LOG("GLRESGenerator", pWARN) << "p = [ " << p4.Px() << " , " << p4.Py() << " , " << p4.Pz() << " ]";
235  LOG("GLRESGenerator", pWARN) << "m = " << p4.M() << " // mpdg = " << massPDG;
236  p4.SetXYZT(0,0,0,massPDG);
237  }
238 
239  // copy final state particles to the event record
241 
242  // fix numbering scheme used for mother/daughter assignments
243  int firstmother = -1;
244  int lastmother = -1;
245  int firstchild = -1;
246  int lastchild = -1;
247 
248  if ( parent < 10 ) { // I=10 is the position were W boson is stored
249  if ( TMath::Abs(pdgc)<7 ) { //outgoing quarks: mother will be the boson (saved in position 4)
250  firstmother = 4;
251  firstchild = p->GetFirstChild() - 6; //boson is stored in I=10 with pythia and I=4 for GENIE => shift 6
252  lastchild = p->GetLastChild() - 6;
253  }
254  else if ( TMath::Abs(pdgc)==24 ) { //produced W boson: mother will be the incoming neutrino
255  firstmother = 0;
256  firstchild = p->GetFirstChild() - 6;
257  lastchild = p->GetLastChild() - 6;
258  }
259  else if ( pdgc==22 ) { //radiative photons: mother will be the incoming electron
260  firstmother = 2;
261  }
262  }
263  else { //rest
264  firstmother = parent - 6; //shift to match boson position
265  firstchild = (p->GetFirstChild()==0) ? p->GetFirstChild() - 1 : p->GetFirstChild() - 6;
266  lastchild = (p->GetLastChild()==0) ? p->GetLastChild() - 1 : p->GetLastChild() - 6;
267  }
268 
269  double vx = nu->X4()->X() + p->GetVx()*1e12; //pythia gives position in [mm] while genie uses [fm]
270  double vy = nu->X4()->Y() + p->GetVy()*1e12;
271  double vz = nu->X4()->Z() + p->GetVz()*1e12;
272  double vt = nu->X4()->T() + p->GetTime()*(units::millimeter/units::second);
273  TLorentzVector pos( vx, vy, vz, vt );
274 
275  event->AddParticle(pdgc, ist, firstmother, lastmother, firstchild, lastchild, p4, pos );
276 
277  }
278 
279  }
280 
281 #else
282  LOG("GLRESGenerator", pFATAL)
283  << "Calling GENIE/PYTHIA6 without enabling PYTHIA6";
284  gAbortingInErr = true;
285  std::exit(1);
286 #endif
287 
288 }
double beta(double KE, const simb::MCParticle *part)
TRandom3 & RndLep(void) const
rnd number generator used by final state primary lepton generators
Definition: RandomGen.h:62
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1064
virtual GHepParticle * HitElectron(void) const
Definition: GHepRecord.cxx:316
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
const TLorentzVector * P4(void) const
Definition: GHepParticle.h:78
const int kPdgAntiNuE
Definition: PDGCodes.h:29
const int kPdgWM
Definition: PDGCodes.h:192
#define pFATAL
Definition: Messenger.h:56
static const double kElectronMass
Definition: Constants.h:70
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
int FSPrimLeptonPdg(void) const
final state primary lepton pdg
virtual GHepParticle * Probe(void) const
Definition: GHepRecord.cxx:277
double y(bool selected=false) const
Definition: Kinematics.cxx:112
int Pdg(void) const
Definition: GHepParticle.h:63
Summary information for an interaction.
Definition: Interaction.h:56
static constexpr double millimeter
Definition: Units.h:41
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
TPythia6 * fPythia
PYTHIA6 wrapper class.
static constexpr double second
Definition: Units.h:37
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
const Kinematics & Kine(void) const
Definition: Interaction.h:71
bool IsTau(int pdgc)
Definition: PDGUtils.cxx:205
p
Definition: test.py:223
#define pINFO
Definition: Messenger.h:62
bool IsMuon(int pdgc)
Definition: PDGUtils.cxx:195
const int kPdgAntiNuTau
Definition: PDGCodes.h:33
const int kPdgAntiNuMu
Definition: PDGCodes.h:31
#define pWARN
Definition: Messenger.h:60
static const double kMw
Definition: Constants.h:91
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
double massPDG(int PDGcode)
Definition: includeROOT.h:203
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
const TLorentzVector * X4(void) const
Definition: GHepParticle.h:79
const InitialState & InitState(void) const
Definition: Interaction.h:69
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
bool gAbortingInErr
Definition: Messenger.cxx:34
double ProbeE(RefFrame_t rf) const
static const double kPi
Definition: Constants.h:37
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
def parent(G, child, parent_type)
Definition: graph.py:67
bool IsElectron(int pdgc)
Definition: PDGUtils.cxx:185
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
enum genie::EGHepStatus GHepStatus_t
Initial State information.
Definition: InitialState.h:48

Member Data Documentation

TPythia6* genie::GLRESGenerator::fPythia
mutableprivate

PYTHIA6 wrapper class.

Definition at line 51 of file GLRESGenerator.h.

double genie::GLRESGenerator::fWmin
private

Definition at line 54 of file GLRESGenerator.h.


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