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