gparticle
Loading...
Searching...
No Matches
Gparticle Class Reference

Lightweight particle specification and primary vertex shooter. More...

#include <gparticle.h>

Public Member Functions

 Gparticle (const std::string &name, int multiplicity, double p, double delta_p, const std::string &randomMomentumModel, double theta, double delta_theta, const std::string &thetaModel, double phi, double delta_phi, double avx, double avy, double avz, double adelta_vx, double adelta_vy, double adelta_vz, const std::string &randomVertexModel, const std::shared_ptr< GLogger > &logger, int generator_type=1)
 Constructs a particle configuration from pre-converted G4-unit values.
 
 Gparticle (const Gparticle &)=delete
 Copying is disabled.
 
Gparticleoperator= (const Gparticle &)=delete
 Copy assignment is disabled.
 
 ~Gparticle ()
 Destructor emits a debug-level lifecycle message.
 
void shootParticle (G4ParticleGun *particleGun, G4Event *anEvent)
 Shoots this particle configuration into a Geant4 event.
 
const std::string & getName () const
 Returns the particle name stored in this generator definition.
 
int getPid () const
 Returns the resolved PDG particle id.
 
int getMultiplicity () const
 Returns how many copies are generated per event.
 
double getMomentum () const
 Returns the nominal momentum magnitude.
 
double getTheta () const
 Returns the nominal polar angle.
 
double getPhi () const
 Returns the nominal azimuthal angle.
 
const G4ThreeVector & getVertex () const
 Returns the nominal vertex position.
 
double getDeltaMomentum () const
 Returns the momentum spread parameter (GEMC internal units, MeV).
 
std::string getMomentumModel () const
 Returns the momentum randomization model as a string token.
 
double getDeltaTheta () const
 Returns the polar-angle spread parameter (GEMC internal units, radians).
 
std::string getThetaModel () const
 Returns the theta randomization model as a string token.
 
double getDeltaPhi () const
 Returns the azimuthal-angle spread parameter (GEMC internal units, radians).
 
const G4ThreeVector & getDeltaVertex () const
 Returns the vertex spread vector (GEMC internal units, mm).
 
std::string getVertexModel () const
 Returns the vertex randomization model as a string token.
 
int getGeneratorType () const
 Returns the generator source type.
 
const std::vector< GparticleRuntimeRecord > & getRuntimeRecords () const
 Returns the runtime records from the most recent shoot.
 
void setName (const std::string &n)
 
void setMultiplicity (int m)
 
void setMomentum (double p_mev)
 
void setDeltaMomentum (double dp)
 
void setMomentumModel (const std::string &s)
 
void setTheta (double t_rad)
 
void setDeltaTheta (double dt)
 
void setThetaModel (const std::string &s)
 
void setPhi (double ph_rad)
 
void setDeltaPhi (double dp)
 
void setVertexX (double x)
 
void setVertexY (double y)
 
void setVertexZ (double z)
 
void setDeltaVertexX (double x)
 
void setDeltaVertexY (double y)
 
void setDeltaVertexZ (double z)
 
void setVertexModel (const std::string &s)
 

Static Public Member Functions

static std::shared_ptr< Gparticlecreate_default_gparticle (const std::shared_ptr< GLogger > &log)
 Creates a minimal default particle configuration.
 

Friends

std::ostream & operator<< (std::ostream &os, const Gparticle &gp)
 Stream insertion operator used for pretty-printing configuration summaries.
 
std::ostream & operator<< (std::ostream &os, const std::shared_ptr< Gparticle > &ptr)
 Stream insertion operator overload for shared pointers.
 

Detailed Description

A Gparticle instance represents a generator-level particle configuration that can be used to produce primary vertices in a G4Event through a G4ParticleGun.

The class stores:

  • Identity: particle name and resolved PDG id
  • Multiplicity: number of copies shot per event
  • Kinematics: momentum magnitude and angular parameters
  • Vertex: position and optional spread/randomization
  • Randomization models: selection of uniform/gaussian/cosine (angles) and sphere (vertex)

Configuration is typically created by option parsing utilities and then used during event generation by calling shootParticle().

Logging:

  • A logger is provided at construction and retained for diagnostics.
  • Verbosity 2 typically prints a full configuration summary via the stream operator.

Definition at line 88 of file gparticle.h.

Constructor & Destructor Documentation

◆ Gparticle() [1/2]

Gparticle::Gparticle ( const std::string & name,
int multiplicity,
double p,
double delta_p,
const std::string & randomMomentumModel,
double theta,
double delta_theta,
const std::string & thetaModel,
double phi,
double delta_phi,
double avx,
double avy,
double avz,
double adelta_vx,
double adelta_vy,
double adelta_vz,
const std::string & randomVertexModel,
const std::shared_ptr< GLogger > & logger,
int generator_type = 1 )

All numeric parameters must already be in Geant4 internal units (MeV for momentum, radians for angles, mm for lengths) as returned by gutilities::getG4Number(). Callers — typically option parsers or file readers — are responsible for the conversion.

The particle PDG id is resolved at construction time by consulting the G4ParticleTable using the provided particle name.

Parameters
nameParticle name as understood by G4ParticleTable (e.g. "e-").
multiplicityNumber of particles to shoot per event.
pNominal momentum magnitude in G4 internal units (MeV).
delta_pSpread parameter for momentum randomization (same units as p).
randomMomentumModelRandom model name for momentum (e.g. "uniform", "gaussian").
thetaNominal polar angle in G4 internal units (radians).
delta_thetaSpread parameter for theta (same units as theta).
thetaModelRandom model name for theta (e.g. "uniform", "cosine", "gaussian").
phiNominal azimuthal angle in G4 internal units (radians).
delta_phiSpread parameter for phi (same units as phi).
avxNominal vertex x component in G4 internal units (mm).
avyNominal vertex y component in G4 internal units (mm).
avzNominal vertex z component in G4 internal units (mm).
adelta_vxSpread parameter for vertex x (same units as avx).
adelta_vySpread parameter for vertex y (same units as avy).
adelta_vzSpread parameter for vertex z (same units as avz).
randomVertexModelRandom model name for vertex (e.g. "uniform", "gaussian", "sphere").
loggerLogger used for diagnostics and error reporting.
generator_typeGenerator source type associated with this particle. Inline -gparticle definitions use the default value 1. File-backed particles preserve the source-file type field when the format provides one, such as the Lund type column.

Definition at line 17 of file gparticle.cc.

◆ Gparticle() [2/2]

Gparticle::Gparticle ( const Gparticle & )
delete

The class is intended to be managed through std::shared_ptr and not copied.

◆ ~Gparticle()

Gparticle::~Gparticle ( )
inline

Definition at line 160 of file gparticle.h.

Member Function Documentation

◆ create_default_gparticle()

static std::shared_ptr< Gparticle > Gparticle::create_default_gparticle ( const std::shared_ptr< GLogger > & log)
inlinestatic

This helper returns an electron with:

  • momentum 1 GeV, no spread
  • angles 0 deg, no spread
  • vertex at (0,0,0) cm, no spread
  • uniform random models where relevant
Parameters
logLogger to associate with the constructed particle.
Returns
A shared pointer to the default Gparticle configuration.

Definition at line 194 of file gparticle.h.

◆ getDeltaMomentum()

double Gparticle::getDeltaMomentum ( ) const
inline

Definition at line 267 of file gparticle.h.

◆ getDeltaPhi()

double Gparticle::getDeltaPhi ( ) const
inline

Definition at line 283 of file gparticle.h.

◆ getDeltaTheta()

double Gparticle::getDeltaTheta ( ) const
inline

Definition at line 275 of file gparticle.h.

◆ getDeltaVertex()

const G4ThreeVector & Gparticle::getDeltaVertex ( ) const
inline

Definition at line 286 of file gparticle.h.

◆ getGeneratorType()

int Gparticle::getGeneratorType ( ) const
inline

For inline -gparticle entries this is normally 1. For file-backed entries this preserves the source format's type field; in Lund files, type 1 is the subset propagated to Geant4 and included in the generated_tracked output bank.

Returns
Generator source type.

Definition at line 303 of file gparticle.h.

◆ getMomentum()

double Gparticle::getMomentum ( ) const
inline
Returns
Momentum after unit conversion to GEMC internal units.

Definition at line 243 of file gparticle.h.

◆ getMomentumModel()

std::string Gparticle::getMomentumModel ( ) const
inline

Definition at line 270 of file gparticle.h.

◆ getMultiplicity()

int Gparticle::getMultiplicity ( ) const
inline
Returns
Configured particle multiplicity.

Definition at line 236 of file gparticle.h.

◆ getName()

const std::string & Gparticle::getName ( ) const
inline
Returns
Particle name used for Geant4 lookup and generated-particle output.

Definition at line 222 of file gparticle.h.

◆ getPhi()

double Gparticle::getPhi ( ) const
inline
Returns
Phi after unit conversion to GEMC internal angular units.

Definition at line 257 of file gparticle.h.

◆ getPid()

int Gparticle::getPid ( ) const
inline
Returns
PDG id resolved from G4ParticleTable during construction.

Definition at line 229 of file gparticle.h.

◆ getRuntimeRecords()

const std::vector< GparticleRuntimeRecord > & Gparticle::getRuntimeRecords ( ) const
inline

One record is stored for each generated primary, so a particle with multiplicity N contributes N runtime records.

Definition at line 311 of file gparticle.h.

◆ getTheta()

double Gparticle::getTheta ( ) const
inline
Returns
Theta after unit conversion to GEMC internal angular units.

Definition at line 250 of file gparticle.h.

◆ getThetaModel()

std::string Gparticle::getThetaModel ( ) const
inline

Definition at line 278 of file gparticle.h.

◆ getVertex()

const G4ThreeVector & Gparticle::getVertex ( ) const
inline
Returns
Vertex after unit conversion to GEMC internal length units.

Definition at line 264 of file gparticle.h.

◆ getVertexModel()

std::string Gparticle::getVertexModel ( ) const
inline

Definition at line 289 of file gparticle.h.

◆ operator=()

Gparticle & Gparticle::operator= ( const Gparticle & )
delete

◆ setDeltaMomentum()

void Gparticle::setDeltaMomentum ( double dp)
inline

Definition at line 319 of file gparticle.h.

◆ setDeltaPhi()

void Gparticle::setDeltaPhi ( double dp)
inline

Definition at line 325 of file gparticle.h.

◆ setDeltaTheta()

void Gparticle::setDeltaTheta ( double dt)
inline

Definition at line 322 of file gparticle.h.

◆ setDeltaVertexX()

void Gparticle::setDeltaVertexX ( double x)
inline

Definition at line 329 of file gparticle.h.

◆ setDeltaVertexY()

void Gparticle::setDeltaVertexY ( double y)
inline

Definition at line 330 of file gparticle.h.

◆ setDeltaVertexZ()

void Gparticle::setDeltaVertexZ ( double z)
inline

Definition at line 331 of file gparticle.h.

◆ setMomentum()

void Gparticle::setMomentum ( double p_mev)
inline

Definition at line 318 of file gparticle.h.

◆ setMomentumModel()

void Gparticle::setMomentumModel ( const std::string & s)
inline

Definition at line 320 of file gparticle.h.

◆ setMultiplicity()

void Gparticle::setMultiplicity ( int m)
inline

Definition at line 317 of file gparticle.h.

◆ setName()

void Gparticle::setName ( const std::string & n)
inline

Definition at line 316 of file gparticle.h.

◆ setPhi()

void Gparticle::setPhi ( double ph_rad)
inline

Definition at line 324 of file gparticle.h.

◆ setTheta()

void Gparticle::setTheta ( double t_rad)
inline

Definition at line 321 of file gparticle.h.

◆ setThetaModel()

void Gparticle::setThetaModel ( const std::string & s)
inline

Definition at line 323 of file gparticle.h.

◆ setVertexModel()

void Gparticle::setVertexModel ( const std::string & s)
inline

Definition at line 332 of file gparticle.h.

◆ setVertexX()

void Gparticle::setVertexX ( double x)
inline

Definition at line 326 of file gparticle.h.

◆ setVertexY()

void Gparticle::setVertexY ( double y)
inline

Definition at line 327 of file gparticle.h.

◆ setVertexZ()

void Gparticle::setVertexZ ( double z)
inline

Definition at line 328 of file gparticle.h.

◆ shootParticle()

void Gparticle::shootParticle ( G4ParticleGun * particleGun,
G4Event * anEvent )

The method resolves the particle definition from the G4ParticleTable (using the stored particle name) and then, for each copy defined by the configured multiplicity:

  • sets the particle kinetic energy based on randomized momentum and mass
  • sets the momentum direction based on randomized theta/phi
  • sets the vertex position based on the configured vertex model
  • calls GeneratePrimaryVertex on the provided G4ParticleGun

Error handling:

  • If the G4ParticleTable is unavailable, an error is logged.
  • If the particle is not found, an error is logged.
Parameters
particleGunThe G4ParticleGun used to generate the primary vertices.
anEventThe G4Event that receives the generated primary vertices.

Definition at line 60 of file gparticle.cc.

Friends And Related Symbol Documentation

◆ operator<< [1/2]

std::ostream & operator<< ( std::ostream & os,
const Gparticle & gp )
friend

Definition at line 230 of file gparticle.cc.

◆ operator<< [2/2]

std::ostream & operator<< ( std::ostream & os,
const std::shared_ptr< Gparticle > & ptr )
friend
Parameters
osOutput stream.
ptrPointer to a Gparticle; prints "<null Gparticle>" if null.
Returns
Output stream.

Definition at line 399 of file gparticle.h.


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