gparticle
Loading...
Searching...
No Matches
gparticle.h
Go to the documentation of this file.
1#pragma once
2
27// gemc
28#include <gemc/glogging/glogger.h>
29#include <gemc/guts/gutilities.h>
30
31// geant4
32#include "G4ThreeVector.hh"
33#include "G4ParticleGun.hh"
34
35// c++
36#include <string>
37#include <vector>
38
46{
47 std::string name;
48 int pid = 0;
49 int type = 1;
50 double p = 0;
51 double theta = 0;
52 double phi = 0;
53 G4ThreeVector vertex;
54};
55
56
89{
90public:
125 Gparticle(const std::string& name,
126 int multiplicity,
127 double p,
128 double delta_p,
129 const std::string& randomMomentumModel,
130 double theta,
131 double delta_theta,
132 const std::string& thetaModel,
133 double phi,
134 double delta_phi,
135 double avx,
136 double avy,
137 double avz,
138 double adelta_vx,
139 double adelta_vy,
140 double adelta_vz,
141 const std::string& randomVertexModel,
142 const std::shared_ptr<GLogger>& logger,
143 int generator_type = 1);
144
150 Gparticle(const Gparticle&) = delete;
151
155 Gparticle& operator=(const Gparticle&) = delete;
156
160 ~Gparticle() { log->debug(DESTRUCTOR, "Gparticle"); }
161
180 void shootParticle(G4ParticleGun* particleGun, G4Event* anEvent);
181
194 static std::shared_ptr<Gparticle> create_default_gparticle(const std::shared_ptr<GLogger>& log) {
195 return std::make_shared<Gparticle>(
196 "e-",
197 1,
199 0,
200 "uniform",
201 0,
202 0,
203 "uniform",
204 0,
205 0,
206 0,
207 0,
208 0,
209 0,
210 0,
211 0,
212 "uniform",
213 log
214 );
215 }
216
222 [[nodiscard]] const std::string& getName() const { return name; }
223
229 [[nodiscard]] int getPid() const { return pid; }
230
236 [[nodiscard]] int getMultiplicity() const { return multiplicity; }
237
243 [[nodiscard]] double getMomentum() const { return p; }
244
250 [[nodiscard]] double getTheta() const { return theta; }
251
257 [[nodiscard]] double getPhi() const { return phi; }
258
264 [[nodiscard]] const G4ThreeVector& getVertex() const { return v; }
265
267 [[nodiscard]] double getDeltaMomentum() const { return delta_p; }
268
270 [[nodiscard]] std::string getMomentumModel() const {
271 return gutilities::to_string(randomMomentumModel);
272 }
273
275 [[nodiscard]] double getDeltaTheta() const { return delta_theta; }
276
278 [[nodiscard]] std::string getThetaModel() const {
279 return gutilities::to_string(randomThetaModel);
280 }
281
283 [[nodiscard]] double getDeltaPhi() const { return delta_phi; }
284
286 [[nodiscard]] const G4ThreeVector& getDeltaVertex() const { return delta_v; }
287
289 [[nodiscard]] std::string getVertexModel() const {
290 return gutilities::to_string(randomVertexModel);
291 }
292
303 [[nodiscard]] int getGeneratorType() const { return generator_type; }
304
311 [[nodiscard]] const std::vector<GparticleRuntimeRecord>& getRuntimeRecords() const {
312 return runtimeRecords;
313 }
314
315 // ---- Setters (used by the GUI to reflect user edits in real time) ----
316 void setName(const std::string& n) { name = n; pid = get_pdg_id(); }
317 void setMultiplicity(int m) { multiplicity = m; }
318 void setMomentum(double p_mev) { p = p_mev; }
319 void setDeltaMomentum(double dp) { delta_p = dp; }
320 void setMomentumModel(const std::string& s) { randomMomentumModel = gutilities::stringToRandomModel(s); }
321 void setTheta(double t_rad) { theta = t_rad; }
322 void setDeltaTheta(double dt) { delta_theta = dt; }
323 void setThetaModel(const std::string& s) { randomThetaModel = gutilities::stringToRandomModel(s); }
324 void setPhi(double ph_rad) { phi = ph_rad; }
325 void setDeltaPhi(double dp) { delta_phi = dp; }
326 void setVertexX(double x) { v.setX(x); }
327 void setVertexY(double y) { v.setY(y); }
328 void setVertexZ(double z) { v.setZ(z); }
329 void setDeltaVertexX(double x) { delta_v.setX(x); }
330 void setDeltaVertexY(double y) { delta_v.setY(y); }
331 void setDeltaVertexZ(double z) { delta_v.setZ(z); }
332 void setVertexModel(const std::string& s) { randomVertexModel = gutilities::stringToRandomModel(s); }
333
334private:
336 std::string name;
337
339 int generator_type;
340
342 int pid;
343
345 int multiplicity;
346
348 double p;
349
351 double delta_p;
352
354 gutilities::randomModel randomMomentumModel;
355
357 double theta;
358
360 double delta_theta;
361
363 gutilities::randomModel randomThetaModel;
364
366 double phi;
367
369 double delta_phi;
370
372 G4ThreeVector v;
373
375 G4ThreeVector delta_v;
376
378 gutilities::randomModel randomVertexModel;
379
381 std::shared_ptr<GLogger> log;
382
384 std::vector<GparticleRuntimeRecord> runtimeRecords;
385
387 [[nodiscard]] double get_mass() const;
388
390 friend std::ostream& operator<<(std::ostream& os, const Gparticle& gp);
391
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>";
402 }
403
412 [[nodiscard]] double randomizeNumberFromSigmaWithModel(double center,
413 double delta,
414 gutilities::randomModel model) const;
415
417 double calculateMomentum();
418
424 double calculateKinEnergy(double mass);
425
427 G4ThreeVector calculateBeamDirection(double thetaRad, double phiRad);
428
430 G4ThreeVector calculateVertex();
431
433 int get_pdg_id();
434
435
436 // assigned at shootParticle for each multiplicity
437 double kinenergy;
438 G4ThreeVector beamDirection;
439 G4ThreeVector vertex;
440
441 void setRunTimeQuantities(double ke, G4ThreeVector bd, G4ThreeVector v) {
442 kinenergy = ke;
443 beamDirection = bd;
444 vertex = v;
445 }
446};
447
453using GparticlePtr = std::shared_ptr<Gparticle>;
void debug(debug_type type, Args &&... args) const
Lightweight particle specification and primary vertex shooter.
Definition gparticle.h:89
double getDeltaPhi() const
Returns the azimuthal-angle spread parameter (GEMC internal units, radians).
Definition gparticle.h:283
const G4ThreeVector & getDeltaVertex() const
Returns the vertex spread vector (GEMC internal units, mm).
Definition gparticle.h:286
void setDeltaMomentum(double dp)
Definition gparticle.h:319
void shootParticle(G4ParticleGun *particleGun, G4Event *anEvent)
Shoots this particle configuration into a Geant4 event.
Definition gparticle.cc:60
void setVertexX(double x)
Definition gparticle.h:326
const std::string & getName() const
Returns the particle name stored in this generator definition.
Definition gparticle.h:222
std::string getMomentumModel() const
Returns the momentum randomization model as a string token.
Definition gparticle.h:270
int getPid() const
Returns the resolved PDG particle id.
Definition gparticle.h:229
int getMultiplicity() const
Returns how many copies are generated per event.
Definition gparticle.h:236
friend std::ostream & operator<<(std::ostream &os, const Gparticle &gp)
Stream insertion operator used for pretty-printing configuration summaries.
Definition gparticle.cc:230
void setVertexY(double y)
Definition gparticle.h:327
double getPhi() const
Returns the nominal azimuthal angle.
Definition gparticle.h:257
Gparticle(const Gparticle &)=delete
Copying is disabled.
void setName(const std::string &n)
Definition gparticle.h:316
double getMomentum() const
Returns the nominal momentum magnitude.
Definition gparticle.h:243
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
void setDeltaTheta(double dt)
Definition gparticle.h:322
int getGeneratorType() const
Returns the generator source type.
Definition gparticle.h:303
void setPhi(double ph_rad)
Definition gparticle.h:324
void setMomentum(double p_mev)
Definition gparticle.h:318
double getDeltaTheta() const
Returns the polar-angle spread parameter (GEMC internal units, radians).
Definition gparticle.h:275
void setDeltaVertexX(double x)
Definition gparticle.h:329
void setTheta(double t_rad)
Definition gparticle.h:321
void setVertexZ(double z)
Definition gparticle.h:328
void setDeltaVertexY(double y)
Definition gparticle.h:330
double getTheta() const
Returns the nominal polar angle.
Definition gparticle.h:250
const std::vector< GparticleRuntimeRecord > & getRuntimeRecords() const
Returns the runtime records from the most recent shoot.
Definition gparticle.h:311
void setThetaModel(const std::string &s)
Definition gparticle.h:323
friend std::ostream & operator<<(std::ostream &os, const std::shared_ptr< Gparticle > &ptr)
Stream insertion operator overload for shared pointers.
Definition gparticle.h:399
static std::shared_ptr< Gparticle > create_default_gparticle(const std::shared_ptr< GLogger > &log)
Creates a minimal default particle configuration.
Definition gparticle.h:194
~Gparticle()
Destructor emits a debug-level lifecycle message.
Definition gparticle.h:160
void setMomentumModel(const std::string &s)
Definition gparticle.h:320
const G4ThreeVector & getVertex() const
Returns the nominal vertex position.
Definition gparticle.h:264
double getDeltaMomentum() const
Returns the momentum spread parameter (GEMC internal units, MeV).
Definition gparticle.h:267
void setVertexModel(const std::string &s)
Definition gparticle.h:332
Gparticle & operator=(const Gparticle &)=delete
Copy assignment is disabled.
void setMultiplicity(int m)
Definition gparticle.h:317
std::string getVertexModel() const
Returns the vertex randomization model as a string token.
Definition gparticle.h:289
void setDeltaPhi(double dp)
Definition gparticle.h:325
void setDeltaVertexZ(double z)
Definition gparticle.h:331
std::string getThetaModel() const
Returns the theta randomization model as a string token.
Definition gparticle.h:278
DESTRUCTOR
std::shared_ptr< Gparticle > GparticlePtr
Shared pointer type used for Gparticle instances.
Definition gparticle.h:453
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.
Definition gparticle.h:46
G4ThreeVector vertex
Definition gparticle.h:53