G4StationPrimaryGenerator.cc
Go to the documentation of this file.
2 #include "G4StationSimulator.h"
3 #include "G4Utils.h"
4 #include <utl/Particle.h>
5 #include <utl/AugerUnits.h>
6 
7 #include <G4Event.hh>
8 #include <G4ParticleGun.hh>
9 #include <G4ParticleDefinition.hh>
10 #include <G4ParticleTable.hh>
11 #include <G4SystemOfUnits.hh>
12 
13 #include <sstream>
14 
15 
16 namespace G4StationSimulatorOG {
17 
19  fParticleGun(new G4ParticleGun(1)),
20  fParticleTable(G4ParticleTable::GetParticleTable())
21  { }
22 
23 
25  {
26  delete fParticleGun;
27  }
28 
29 
30  void
32  {
33  const utl::Particle& currentParticle = G4StationSimulator::fgCurrent.GetParticle();
34  G4ParticleDefinition* const particleDef = fParticleTable->FindParticle(currentParticle.GetName());
35  if (!particleDef) {
36  std::ostringstream msg;
37  msg << "Undefined particle type: " << currentParticle.GetName();
38  WARNING(msg);
39  return;
40  }
41 
42  const sdet::Station& currentDetectorStation = G4StationSimulator::fgCurrent.GetDetectorStation();
43  const auto stationCS = currentDetectorStation.GetLocalCoordinateSystem();
44  const G4ThreeVector position = ToG4Vector(currentParticle.GetPosition(), stationCS, meter/utl::meter);
45  const G4ThreeVector direction = ToG4Vector(currentParticle.GetDirection(), stationCS, 1);
46 
47  // set properties and pass to particle gun
48  fParticleGun->SetParticleDefinition(particleDef);
49  fParticleGun->SetParticlePosition(position);
50  fParticleGun->SetParticleMomentumDirection(direction);
51 
52  // Convert from Auger units to G4 units
53  fParticleGun->SetParticleEnergy(currentParticle.GetKineticEnergy() * (CLHEP::eV / utl::eV));
54  fParticleGun->SetParticleTime(currentParticle.GetTime().GetInterval());
55  fParticleGun->GeneratePrimaryVertex(event);
56  }
57 
58 }
const double eV
Definition: GalacticUnits.h:35
constexpr double eV
Definition: AugerUnits.h:185
const TimeInterval & GetTime() const
Arrival time delay at the ground.
Definition: Particle.h:122
Detector description interface for Station-related data.
std::string GetName() const
string with particle name
Definition: Particle.h:104
Describes a particle for Simulation.
Definition: Particle.h:26
const double meter
Definition: GalacticUnits.h:29
constexpr double meter
Definition: AugerUnits.h:81
utl::CoordinateSystemPtr GetLocalCoordinateSystem() const
Get the Auger reference system centered on the tank.
#define WARNING(message)
Macro for logging warning messages.
Definition: ErrorLogger.h:163
double GetInterval() const
Get the time interval as a double (in Auger base units)
Definition: TimeInterval.h:69
double GetKineticEnergy() const
Get kinetic energy of the particle.
Definition: Particle.h:130
G4ThreeVector ToG4Vector(const V &v, const utl::CoordinateSystemPtr &cs, const double unitConversion)
Definition: G4Utils.h:15
const Point & GetPosition() const
Position of the particle.
Definition: Particle.h:110
const Vector & GetDirection() const
Unit vector giving particle direction.
Definition: Particle.h:114

, generated on Tue Sep 26 2023.