28#include <gemc/glogging/glogger.h>
29#include <gemc/guts/gutilities.h>
32#include "G4ThreeVector.hh"
33#include "G4ParticleGun.hh"
129 const std::string& randomMomentumModel,
132 const std::string& thetaModel,
141 const std::string& randomVertexModel,
142 const std::shared_ptr<GLogger>& logger,
143 int generator_type = 1);
180 void shootParticle(G4ParticleGun* particleGun, G4Event* anEvent);
195 return std::make_shared<Gparticle>(
222 [[nodiscard]]
const std::string&
getName()
const {
return name; }
229 [[nodiscard]]
int getPid()
const {
return pid; }
250 [[nodiscard]]
double getTheta()
const {
return theta; }
257 [[nodiscard]]
double getPhi()
const {
return phi; }
264 [[nodiscard]]
const G4ThreeVector&
getVertex()
const {
return v; }
312 return runtimeRecords;
316 void setName(
const std::string& n) { name = n; pid = get_pdg_id(); }
324 void setPhi(
double ph_rad) { phi = ph_rad; }
375 G4ThreeVector delta_v;
381 std::shared_ptr<GLogger> log;
384 std::vector<GparticleRuntimeRecord> runtimeRecords;
387 [[nodiscard]]
double get_mass()
const;
399 friend std::ostream&
operator<<(std::ostream& os,
const std::shared_ptr<Gparticle>& ptr) {
400 if (ptr)
return os << *ptr;
401 else return os <<
"<null Gparticle>";
412 [[nodiscard]]
double randomizeNumberFromSigmaWithModel(
double center,
417 double calculateMomentum();
424 double calculateKinEnergy(
double mass);
427 G4ThreeVector calculateBeamDirection(
double thetaRad,
double phiRad);
430 G4ThreeVector calculateVertex();
438 G4ThreeVector beamDirection;
439 G4ThreeVector vertex;
441 void setRunTimeQuantities(
double ke, G4ThreeVector bd, G4ThreeVector v) {
void debug(debug_type type, Args &&... args) const
Lightweight particle specification and primary vertex shooter.
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).
void setDeltaMomentum(double dp)
void shootParticle(G4ParticleGun *particleGun, G4Event *anEvent)
Shoots this particle configuration into a Geant4 event.
void setVertexX(double x)
const std::string & getName() const
Returns the particle name stored in this generator definition.
std::string getMomentumModel() const
Returns the momentum randomization model as a string token.
int getPid() const
Returns the resolved PDG particle id.
int getMultiplicity() const
Returns how many copies are generated per event.
friend std::ostream & operator<<(std::ostream &os, const Gparticle &gp)
Stream insertion operator used for pretty-printing configuration summaries.
void setVertexY(double y)
double getPhi() const
Returns the nominal azimuthal angle.
Gparticle(const Gparticle &)=delete
Copying is disabled.
void setName(const std::string &n)
double getMomentum() const
Returns the nominal momentum magnitude.
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.
void setDeltaTheta(double dt)
int getGeneratorType() const
Returns the generator source type.
void setPhi(double ph_rad)
void setMomentum(double p_mev)
double getDeltaTheta() const
Returns the polar-angle spread parameter (GEMC internal units, radians).
void setDeltaVertexX(double x)
void setTheta(double t_rad)
void setVertexZ(double z)
void setDeltaVertexY(double y)
double getTheta() const
Returns the nominal polar angle.
const std::vector< GparticleRuntimeRecord > & getRuntimeRecords() const
Returns the runtime records from the most recent shoot.
void setThetaModel(const std::string &s)
friend std::ostream & operator<<(std::ostream &os, const std::shared_ptr< Gparticle > &ptr)
Stream insertion operator overload for shared pointers.
static std::shared_ptr< Gparticle > create_default_gparticle(const std::shared_ptr< GLogger > &log)
Creates a minimal default particle configuration.
~Gparticle()
Destructor emits a debug-level lifecycle message.
void setMomentumModel(const std::string &s)
const G4ThreeVector & getVertex() const
Returns the nominal vertex position.
double getDeltaMomentum() const
Returns the momentum spread parameter (GEMC internal units, MeV).
void setVertexModel(const std::string &s)
Gparticle & operator=(const Gparticle &)=delete
Copy assignment is disabled.
void setMultiplicity(int m)
std::string getVertexModel() const
Returns the vertex randomization model as a string token.
void setDeltaPhi(double dp)
void setDeltaVertexZ(double z)
std::string getThetaModel() const
Returns the theta randomization model as a string token.
std::shared_ptr< Gparticle > GparticlePtr
Shared pointer type used for Gparticle instances.
double getG4Number(const string &v, bool warnIfNotUnit=false)
randomModel stringToRandomModel(const std::string &str)
constexpr const char * to_string(randomModel m) noexcept
Runtime values for one generated primary particle.