9#include "G4ParticleTable.hh"
10#include "Randomize.hh"
21 const std::string& arandomMomentumModel,
24 const std::string& arandomThetaModel,
33 const std::string& arandomVertexModel,
34 const std::shared_ptr<GLogger>& logger,
35 int agenerator_type) :
37 generator_type(agenerator_type),
38 multiplicity(amultiplicity),
43 delta_theta(adelta_theta),
46 delta_phi(adelta_phi),
47 v(G4ThreeVector(avx, avy, avz)),
48 delta_v(G4ThreeVector(adelta_vx, adelta_vy, adelta_vz)),
61 auto particleTable = G4ParticleTable::GetParticleTable();
62 runtimeRecords.clear();
66 auto particleDef = particleTable->FindParticle(name);
70 double mass = particleDef->GetPDGMass();
71 particleGun->SetParticleDefinition(particleDef);
75 for (
int i = 0; i < multiplicity; i++) {
76 auto pmev = calculateMomentum();
77 auto kenergy = sqrt(pmev * pmev + mass * mass) - mass;
78 auto thetaRad = randomizeNumberFromSigmaWithModel(theta, delta_theta, randomThetaModel) / CLHEP::rad;
79 auto phiRad = randomizeNumberFromSigmaWithModel(phi, delta_phi,
gutilities::uniform) / CLHEP::rad;
80 auto beamDirection = calculateBeamDirection(thetaRad, phiRad);
81 auto vertex = calculateVertex();
83 setRunTimeQuantities(kenergy, beamDirection, vertex);
84 runtimeRecords.push_back({
89 thetaRad * CLHEP::rad,
94 particleGun->SetParticleEnergy(kenergy);
95 particleGun->SetParticleMomentumDirection(beamDirection);
96 particleGun->SetParticlePosition(vertex);
97 particleGun->GeneratePrimaryVertex(anEvent);
104 "Particle <", name,
"> not found in G4ParticleTable* ", particleTable);
109 "G4ParticleTable not found - G4ParticleGun*: ", particleGun);
114double Gparticle::calculateMomentum() {
116 double pmev = randomizeNumberFromSigmaWithModel(p, delta_p, randomMomentumModel);
121double Gparticle::calculateKinEnergy(
double mass) {
122 double pmev = calculateMomentum();
124 return sqrt(pmev * pmev + mass * mass) - mass;
128G4ThreeVector Gparticle::calculateBeamDirection(
double thetaRad,
double phiRad) {
129 G4ThreeVector pdir = G4ThreeVector(
130 cos(phiRad) * sin(thetaRad),
131 sin(phiRad) * sin(thetaRad),
138G4ThreeVector Gparticle::calculateVertex() {
141 switch (randomVertexModel) {
160 double max_radius = delta_v.r();
167 radius = x * x + y * y + z * z;
169 while (radius > max_radius * max_radius);
190double Gparticle::randomizeNumberFromSigmaWithModel(
double center,
double delta,
gutilities::randomModel model)
const {
194 return center + ( 2.0 * G4UniformRand() - 1.0 ) * delta;
198 return G4RandGauss::shoot(center, delta);
204 double lower = ( center - delta ) / CLHEP::rad;
205 double upper = ( center + delta ) / CLHEP::rad;
206 double center_rad = 0;
210 do { center_rad = acos(1 - 2 * G4UniformRand()); }
211 while (center_rad < lower || center_rad > upper);
215 center_rad = theta / CLHEP::rad;
218 return center_rad * CLHEP::rad;
235 constexpr int label_w = 15;
236 constexpr int value_w = 12;
239 auto show = [&](
const std::string& label,
const auto& value) {
240 os << left << setw(label_w) << label <<
' '
241 << setw(value_w) << right << setw(value_w) << value <<
'\n';
245 auto showf = [&](
const std::string& label,
double value,
int prec = 3) {
246 std::streamsize old_prec = os.precision();
247 auto old_flags = os.flags();
249 os << left << setw(label_w) << label <<
' '
250 << right << setw(value_w) << std::fixed
251 << std::setprecision(prec) << value <<
'\n';
253 os.precision(old_prec);
258 auto show_pm = [&](
const std::string& label,
259 double val,
double err,
261 std::streamsize old_prec = os.precision();
262 auto old_flags = os.flags();
264 os << left << setw(label_w) << label <<
' '
265 << right << setw(value_w) << std::fixed << setw(value_w)
266 << std::setprecision(prec) << val
267 <<
" ± " << std::setprecision(prec) << err <<
'\n';
269 os.precision(old_prec);
277 <<
" ┌─────────────────────────────────────────────────┐\n"
278 <<
" │ GParticle │\n"
279 <<
" └─────────────────────────────────────────────────┘\n";
284 os << left << setw(label_w) <<
" name:" << right << setw(value_w)
285 << gp.name <<
"(pid " << gp.pid <<
")\n";
287 show(
" multiplicity:", std::to_string(gp.multiplicity));
288 showf(
" mass [MeV]:", gp.get_mass());
292 gp.delta_p / CLHEP::MeV);
294 show(
" p model:",
to_string(gp.randomMomentumModel));
296 show_pm(
" theta [deg]:",
297 gp.theta / CLHEP::deg,
298 gp.delta_theta / CLHEP::deg);
300 show(
" theta model:",
to_string(gp.randomThetaModel));
302 show_pm(
" phi [deg]:",
304 gp.delta_phi / CLHEP::deg);
306 os << left << setw(label_w) <<
" vertex [cm]:" <<
' '
307 << gp.v <<
" ± " << gp.delta_v <<
'\n';
309 show(
" vertex model:",
to_string(gp.randomVertexModel));
311 showf(
" kinematic energy [MeV]:", gp.kinenergy / CLHEP::MeV);
312 os << left << setw(label_w) <<
" beam direction :" <<
' '
313 << gp.beamDirection <<
'\n';
315 os << left << setw(label_w) <<
" vertex :" <<
' '
316 << gp.vertex <<
'\n';
323int Gparticle::get_pdg_id() {
324 auto particleTable = G4ParticleTable::GetParticleTable();
327 auto particleDef = particleTable->FindParticle(name);
329 if (particleDef !=
nullptr) {
return particleDef->GetPDGEncoding(); }
332 "Particle <", name,
"> not found in G4ParticleTable* ", particleTable);
337 "G4ParticleTable not found - G4ParticleGun*: ", particleTable);
342double Gparticle::get_mass()
const {
343 auto particleTable = G4ParticleTable::GetParticleTable();
346 auto particleDef = particleTable->FindParticle(name);
349 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 &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.
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.
randomModel stringToRandomModel(const std::string &str)
constexpr const char * to_string(randomModel m) noexcept