12#include "G4ParticleTable.hh"
13#include "Randomize.hh"
26 const std::string& punit,
27 const std::string& arandomMomentumModel,
30 const std::string& arandomThetaModel,
33 const std::string& aunit,
40 const std::string& vunit,
41 const std::string& arandomVertexModel,
42 const std::shared_ptr<GLogger>& logger) :
44 multiplicity(amultiplicity),
59 delta_v(G4ThreeVector(
76 auto particleTable = G4ParticleTable::GetParticleTable();
80 auto particleDef = particleTable->FindParticle(name);
84 double mass = particleDef->GetPDGMass();
85 particleGun->SetParticleDefinition(particleDef);
89 for (
int i = 0; i < multiplicity; i++) {
90 auto kenergy = calculateKinEnergy(mass);
91 auto beamDirection = calculateBeamDirection();
92 auto vertex = calculateVertex();
94 setRunTimeQuantities(kenergy, beamDirection, vertex);
96 particleGun->SetParticleEnergy(kenergy);
97 particleGun->SetParticleMomentumDirection(beamDirection);
98 particleGun->SetParticlePosition(vertex);
99 particleGun->GeneratePrimaryVertex(anEvent);
106 "Particle <", name,
"> not found in G4ParticleTable* ", particleTable);
111 "G4ParticleTable not found - G4ParticleGun*: ", particleGun);
116double Gparticle::calculateMomentum() {
118 double pmev = randomizeNumberFromSigmaWithModel(p, delta_p, randomMomentumModel);
123double Gparticle::calculateKinEnergy(
double mass) {
124 double pmev = calculateMomentum();
126 return sqrt(pmev * pmev + mass * mass) - mass;
130G4ThreeVector Gparticle::calculateBeamDirection() {
132 double thetaRad = randomizeNumberFromSigmaWithModel(theta, delta_theta, randomThetaModel) / CLHEP::rad;
133 double phiRad = randomizeNumberFromSigmaWithModel(phi, delta_phi,
gutilities::uniform) / CLHEP::rad;
135 G4ThreeVector pdir = G4ThreeVector(
136 cos(phiRad) * sin(thetaRad),
137 sin(phiRad) * sin(thetaRad),
144G4ThreeVector Gparticle::calculateVertex() {
147 switch (randomVertexModel) {
166 double max_radius = delta_v.r();
173 radius = x * x + y * y + z * z;
175 while (radius > max_radius);
196double Gparticle::randomizeNumberFromSigmaWithModel(
double center,
double delta,
gutilities::randomModel model)
const {
200 return center + ( 2.0 * G4UniformRand() - 1.0 ) * delta;
204 return G4RandGauss::shoot(center, delta);
210 double lower = ( center - delta ) / CLHEP::rad;
211 double upper = ( center + delta ) / CLHEP::rad;
212 double center_rad = 0;
216 do { center_rad = acos(1 - 2 * G4UniformRand()); }
217 while (center_rad < lower || center_rad > upper);
221 center_rad = theta / CLHEP::rad;
224 return center_rad * CLHEP::rad;
241 constexpr int label_w = 15;
242 constexpr int value_w = 12;
245 auto show = [&](
const std::string& label,
const auto& value) {
246 os << left << setw(label_w) << label <<
' '
247 << setw(value_w) << right << setw(value_w) << value <<
'\n';
251 auto showf = [&](
const std::string& label,
double value,
int prec = 3) {
252 std::streamsize old_prec = os.precision();
253 auto old_flags = os.flags();
255 os << left << setw(label_w) << label <<
' '
256 << right << setw(value_w) << std::fixed
257 << std::setprecision(prec) << value <<
'\n';
259 os.precision(old_prec);
264 auto show_pm = [&](
const std::string& label,
265 double val,
double err,
267 std::streamsize old_prec = os.precision();
268 auto old_flags = os.flags();
270 os << left << setw(label_w) << label <<
' '
271 << right << setw(value_w) << std::fixed << setw(value_w)
272 << std::setprecision(prec) << val
273 <<
" ± " << std::setprecision(prec) << err <<
'\n';
275 os.precision(old_prec);
283 <<
" ┌─────────────────────────────────────────────────┐\n"
284 <<
" │ GParticle │\n"
285 <<
" └─────────────────────────────────────────────────┘\n";
290 os << left << setw(label_w) <<
" name:" << right << setw(value_w)
291 << gp.name <<
"(pid " << gp.pid <<
")\n";
293 show(
" multiplicity:", std::to_string(gp.multiplicity));
294 showf(
" mass [MeV]:", gp.get_mass());
298 gp.delta_p / CLHEP::MeV);
300 show(
" p model:", to_string(gp.randomMomentumModel));
302 show_pm(
" theta [deg]:",
303 gp.theta / CLHEP::deg,
304 gp.delta_theta / CLHEP::deg);
306 show(
" theta model:", to_string(gp.randomThetaModel));
308 show_pm(
" phi [deg]:",
310 gp.delta_phi / CLHEP::deg);
312 os << left << setw(label_w) <<
" vertex [cm]:" <<
' '
313 << gp.v <<
" ± " << gp.delta_v <<
'\n';
315 show(
" vertex model:", to_string(gp.randomVertexModel));
317 showf(
" kinematic energy [MeV]:", gp.kinenergy / CLHEP::MeV);
318 os << left << setw(label_w) <<
" beam direction :" <<
' '
319 << gp.beamDirection <<
'\n';
321 os << left << setw(label_w) <<
" vertex :" <<
' '
322 << gp.vertex <<
'\n';
329int Gparticle::get_pdg_id() {
330 auto particleTable = G4ParticleTable::GetParticleTable();
333 auto particleDef = particleTable->FindParticle(name);
335 if (particleDef !=
nullptr) {
return particleDef->GetPDGEncoding(); }
338 "Particle <", name,
"> not found in G4ParticleTable* ", particleTable);
343 "G4ParticleTable not found - G4ParticleGun*: ", particleTable);
348double Gparticle::get_mass()
const {
349 auto particleTable = G4ParticleTable::GetParticleTable();
352 auto particleDef = particleTable->FindParticle(name);
355 double mass = particleDef->GetPDGMass();
void debug(debug_type type, Args &&... args) const
void info(int level, Args &&... args) const
void error(int exit_code, Args &&... args) const
Lightweight particle specification and primary vertex shooter.
void shootParticle(G4ParticleGun *particleGun, G4Event *anEvent)
Shoots this particle configuration into a Geant4 event.
Gparticle(const std::string &name, int multiplicity, double p, double delta_p, const std::string &punit, const std::string &randomMomentumModel, double theta, double delta_theta, const std::string &thetaModel, double phi, double delta_phi, const std::string &aunit, double avx, double avy, double avz, double adelta_vx, double adelta_vy, double adelta_vz, const std::string &vunit, const std::string &randomVertexModel, const std::shared_ptr< GLogger > &logger)
Constructs a particle configuration from user-facing parameters.
Conventions and error codes for the gparticle module.
#define ERR_GPARTICLETABLENOTFOUND
G4ParticleTable could not be obtained (unexpected runtime state).
#define ERR_GPARTICLENOTFOUND
Requested particle name was not found in the G4ParticleTable.
std::ostream & operator<<(std::ostream &os, const Gparticle &gp)
Definition of the Gparticle class used by the gparticle module.
Public API for defining and parsing gparticle-related options.
double getG4Number(const string &v, bool warnIfNotUnit=false)
randomModel stringToRandomModel(const std::string &str)