gparticle
Loading...
Searching...
No Matches
gparticle.cc
Go to the documentation of this file.
1// guts
2#include "gutilities.h"
3
4// gparticle
5#include "gparticle.h"
6#include "gparticle_options.h"
7
8// gemc
10
11// geant4
12#include "G4ParticleTable.hh"
13#include "Randomize.hh"
14
15using std::to_string;
16using std::ostream;
17using std::endl;
18using std::string;
19
20// Constructor based on parameters.
21// Detailed API documentation is in gparticle.h.
22Gparticle::Gparticle(const std::string& aname,
23 int amultiplicity,
24 double ap,
25 double adelta_p,
26 const std::string& punit,
27 const std::string& arandomMomentumModel,
28 double atheta,
29 double adelta_theta,
30 const std::string& arandomThetaModel,
31 double aphi,
32 double adelta_phi,
33 const std::string& aunit,
34 double avx,
35 double avy,
36 double avz,
37 double adelta_vx,
38 double adelta_vy,
39 double adelta_vz,
40 const std::string& vunit,
41 const std::string& arandomVertexModel,
42 const std::shared_ptr<GLogger>& logger) :
43 name(aname),
44 multiplicity(amultiplicity),
45 // Convert user values + unit strings into Geant4 numeric values through gutilities.
46 p(gutilities::getG4Number(to_string(ap) + "*" + punit)),
47 delta_p(gutilities::getG4Number(to_string(adelta_p) + "*" + punit)),
48 randomMomentumModel(gutilities::stringToRandomModel(arandomMomentumModel)),
49 theta(gutilities::getG4Number(to_string(atheta) + "*" + aunit)),
50 delta_theta(gutilities::getG4Number(to_string(adelta_theta) + "*" + aunit)),
51 randomThetaModel(gutilities::stringToRandomModel(arandomThetaModel)),
52 phi(gutilities::getG4Number(to_string(aphi) + "*" + aunit)),
53 delta_phi(gutilities::getG4Number(to_string(adelta_phi) + "*" + aunit)),
54 v(G4ThreeVector(
55 gutilities::getG4Number(to_string(avx) + "*" + vunit),
56 gutilities::getG4Number(to_string(avy) + "*" + vunit),
57 gutilities::getG4Number(to_string(avz) + "*" + vunit)
58 )),
59 delta_v(G4ThreeVector(
60 gutilities::getG4Number(to_string(adelta_vx) + "*" + vunit),
61 gutilities::getG4Number(to_string(adelta_vy) + "*" + vunit),
62 gutilities::getG4Number(to_string(adelta_vz) + "*" + vunit)
63 )),
64 randomVertexModel(gutilities::stringToRandomModel(arandomVertexModel)),
65 log(logger) {
66 // Resolve PDG id immediately so errors are detected early and configuration printing is complete.
67 pid = get_pdg_id();
68
69 log->debug(CONSTRUCTOR, "Gparticle");
70}
71
72
73// Shoots this particle into the provided event using the provided Geant4 particle gun.
74// Detailed API documentation is in gparticle.h.
75void Gparticle::shootParticle(G4ParticleGun* particleGun, G4Event* anEvent) {
76 auto particleTable = G4ParticleTable::GetParticleTable();
77
78 if (particleTable) {
79 // Resolve the particle definition by name.
80 auto particleDef = particleTable->FindParticle(name);
81
82 if (particleDef) {
83 // Mass is used to convert randomized momentum magnitude into kinetic energy.
84 double mass = particleDef->GetPDGMass();
85 particleGun->SetParticleDefinition(particleDef);
86
87
88 // Shoot one primary vertex per multiplicity.
89 for (int i = 0; i < multiplicity; i++) {
90 auto kenergy = calculateKinEnergy(mass);
91 auto beamDirection = calculateBeamDirection();
92 auto vertex = calculateVertex();
93
94 setRunTimeQuantities(kenergy, beamDirection, vertex);
95
96 particleGun->SetParticleEnergy(kenergy);
97 particleGun->SetParticleMomentumDirection(beamDirection);
98 particleGun->SetParticlePosition(vertex);
99 particleGun->GeneratePrimaryVertex(anEvent);
100
101 log->info(2, *this);
102 }
103 }
104 else {
106 "Particle <", name, "> not found in G4ParticleTable* ", particleTable);
107 }
108 }
109 else {
111 "G4ParticleTable not found - G4ParticleGun*: ", particleGun);
112 }
113}
114
115
116double Gparticle::calculateMomentum() {
117 // randomizeNumberFromSigmaWithModel applies the model-dependent interpretation of delta.
118 double pmev = randomizeNumberFromSigmaWithModel(p, delta_p, randomMomentumModel);
119
120 return pmev;
121}
122
123double Gparticle::calculateKinEnergy(double mass) {
124 double pmev = calculateMomentum();
125
126 return sqrt(pmev * pmev + mass * mass) - mass;
127}
128
129
130G4ThreeVector Gparticle::calculateBeamDirection() {
131 // Convert to radians for trigonometric functions.
132 double thetaRad = randomizeNumberFromSigmaWithModel(theta, delta_theta, randomThetaModel) / CLHEP::rad;
133 double phiRad = randomizeNumberFromSigmaWithModel(phi, delta_phi, gutilities::uniform) / CLHEP::rad;
134
135 G4ThreeVector pdir = G4ThreeVector(
136 cos(phiRad) * sin(thetaRad),
137 sin(phiRad) * sin(thetaRad),
138 cos(thetaRad)
139 );
140
141 return pdir;
142}
143
144G4ThreeVector Gparticle::calculateVertex() {
145 double x, y, z;
146
147 switch (randomVertexModel) {
149 // Component-wise uniform sampling around the nominal vertex.
150 x = randomizeNumberFromSigmaWithModel(v.x(), delta_v.x(), gutilities::uniform);
151 y = randomizeNumberFromSigmaWithModel(v.y(), delta_v.y(), gutilities::uniform);
152 z = randomizeNumberFromSigmaWithModel(v.z(), delta_v.z(), gutilities::uniform);
153 break;
154
156 // Component-wise Gaussian sampling around the nominal vertex (deltas used as sigmas).
157 x = randomizeNumberFromSigmaWithModel(v.x(), delta_v.x(), gutilities::gaussian);
158 y = randomizeNumberFromSigmaWithModel(v.y(), delta_v.y(), gutilities::gaussian);
159 z = randomizeNumberFromSigmaWithModel(v.z(), delta_v.z(), gutilities::gaussian);
160 break;
161
162 case gutilities::sphere: {
163 // Sample an offset inside a sphere-like region whose maximum radius is delta_v.r().
164 // Assumes all three components have comparable spread such that delta_v.r() is meaningful.
165 double radius;
166 double max_radius = delta_v.r();
167
168 // Rejection sampling: generate a random point in a cube until it lies within the radius bound.
169 do {
170 x = randomizeNumberFromSigmaWithModel(0, max_radius, gutilities::uniform);
171 y = randomizeNumberFromSigmaWithModel(0, max_radius, gutilities::uniform);
172 z = randomizeNumberFromSigmaWithModel(0, max_radius, gutilities::uniform);
173 radius = x * x + y * y + z * z;
174 }
175 while (radius > max_radius);
176
177 // Offset the sampled point by the nominal vertex.
178 x = x + v.x();
179 y = y + v.y();
180 z = z + v.z();
181 break;
182 }
183
184 default:
185 // Unknown model: fall back to deterministic vertex.
186 x = v.x();
187 y = v.y();
188 z = v.z();
189 break;
190 }
191
192 return {x, y, z};
193}
194
195
196double Gparticle::randomizeNumberFromSigmaWithModel(double center, double delta, gutilities::randomModel model) const {
197 switch (model) {
199 // Uniform in [center-delta, center+delta].
200 return center + ( 2.0 * G4UniformRand() - 1.0 ) * delta;
201
203 // Gaussian with mean=center and sigma=delta.
204 return G4RandGauss::shoot(center, delta);
205
206 case gutilities::cosine: {
207 // assuming this is an angle with corrected units
208 // For cosine-weighted sampling we work in radians, sample theta with sin(theta) weighting,
209 // and then enforce the requested [center-delta, center+delta] range.
210 double lower = ( center - delta ) / CLHEP::rad;
211 double upper = ( center + delta ) / CLHEP::rad;
212 double center_rad = 0;
213
214 if (lower < upper) {
215 // Generate theta such that cos(theta) is uniform, which corresponds to sin(theta) weighting.
216 do { center_rad = acos(1 - 2 * G4UniformRand()); }
217 while (center_rad < lower || center_rad > upper);
218 }
219 else {
220 // Degenerate range: fall back to the stored theta value.
221 center_rad = theta / CLHEP::rad;
222 }
223
224 return center_rad * CLHEP::rad;
225 }
226
227 default:
228 // Unknown model: no randomization.
229 return center;
230 }
231}
232
233// ---------------------------------------------------------------------------
234// pretty printer
235// ---------------------------------------------------------------------------
236std::ostream& operator<<(std::ostream& os, const Gparticle& gp) {
237 using std::left;
238 using std::right;
239 using std::setw;
240
241 constexpr int label_w = 15; // width for the field name (with ':')
242 constexpr int value_w = 12; // width for the main column
243
244 // helper: plain value
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';
248 };
249
250 // helper: double value with N decimals
251 auto showf = [&](const std::string& label, double value, int prec = 3) {
252 std::streamsize old_prec = os.precision(); // save current
253 auto old_flags = os.flags(); // save flags
254
255 os << left << setw(label_w) << label << ' '
256 << right << setw(value_w) << std::fixed
257 << std::setprecision(prec) << value << '\n';
258
259 os.precision(old_prec); // restore
260 os.flags(old_flags);
261 };
262
263 // helper: value ± error (both doubles)
264 auto show_pm = [&](const std::string& label,
265 double val, double err,
266 int prec = 3) {
267 std::streamsize old_prec = os.precision();
268 auto old_flags = os.flags();
269
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';
274
275 os.precision(old_prec);
276 os.flags(old_flags);
277 };
278
279 // -----------------------------------------------------------------------
280 // header block
281 // -----------------------------------------------------------------------
282 os << '\n'
283 << " ┌─────────────────────────────────────────────────┐\n"
284 << " │ GParticle │\n"
285 << " └─────────────────────────────────────────────────┘\n";
286
287 // -----------------------------------------------------------------------
288 // fields
289 // -----------------------------------------------------------------------
290 os << left << setw(label_w) << " name:" << right << setw(value_w)
291 << gp.name << "(pid " << gp.pid << ")\n";
292
293 show(" multiplicity:", std::to_string(gp.multiplicity));
294 showf(" mass [MeV]:", gp.get_mass());
295
296 show_pm(" p [MeV]:",
297 gp.p / CLHEP::MeV,
298 gp.delta_p / CLHEP::MeV);
299
300 show(" p model:", to_string(gp.randomMomentumModel));
301
302 show_pm(" theta [deg]:",
303 gp.theta / CLHEP::deg,
304 gp.delta_theta / CLHEP::deg);
305
306 show(" theta model:", to_string(gp.randomThetaModel));
307
308 show_pm(" phi [deg]:",
309 gp.phi / CLHEP::deg,
310 gp.delta_phi / CLHEP::deg);
311
312 os << left << setw(label_w) << " vertex [cm]:" << ' '
313 << gp.v << " ± " << gp.delta_v << '\n';
314
315 show(" vertex model:", to_string(gp.randomVertexModel));
316
317 showf(" kinematic energy [MeV]:", gp.kinenergy / CLHEP::MeV);
318 os << left << setw(label_w) << " beam direction :" << ' '
319 << gp.beamDirection << '\n';
320
321 os << left << setw(label_w) << " vertex :" << ' '
322 << gp.vertex << '\n';
323
324
325 return os;
326}
327
328
329int Gparticle::get_pdg_id() {
330 auto particleTable = G4ParticleTable::GetParticleTable();
331
332 if (particleTable) {
333 auto particleDef = particleTable->FindParticle(name);
334
335 if (particleDef != nullptr) { return particleDef->GetPDGEncoding(); }
336 else {
338 "Particle <", name, "> not found in G4ParticleTable* ", particleTable);
339 }
340 }
341 else {
343 "G4ParticleTable not found - G4ParticleGun*: ", particleTable);
344 }
345}
346
347
348double Gparticle::get_mass() const {
349 auto particleTable = G4ParticleTable::GetParticleTable();
350
351 if (particleTable) {
352 auto particleDef = particleTable->FindParticle(name);
353
354 if (particleDef) {
355 double mass = particleDef->GetPDGMass();
356 return mass;
357 }
358 }
359 return 0;
360}
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.
Definition gparticle.h:71
void shootParticle(G4ParticleGun *particleGun, G4Event *anEvent)
Shoots this particle configuration into a Geant4 event.
Definition gparticle.cc:75
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.
Definition gparticle.cc:22
CONSTRUCTOR
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 gparticle.cc:236
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)