33#include "G4VPhysicalVolume.hh"
34#include "G4GDMLParser.hh"
35#include "G4Electron.hh"
38#include "G4SystemOfUnits.hh"
42 G4VFastSimulationModel(modelName, envelope) {
49 G4VFastSimulationModel(modelName) {
58 G4VPhysicalVolume* physicalVolume) {
60 G4GDMLParser* parser =
new G4GDMLParser();
61 remove(
"garfieldGeometry.gdml");
62 parser->Write(
"garfieldGeometry.gdml", physicalVolume,
false);
67 const G4ParticleDefinition& particleType) {
68 G4String particleName = particleType.GetParticleName();
76 const G4FastTrack& fastTrack) {
77 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
78 G4String particleName =
79 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
88 G4FastStep& fastStep) {
90 G4TouchableHandle theTouchable =
91 fastTrack.GetPrimaryTrack()->GetTouchableHandle();
92 G4String name = theTouchable->GetVolume()->GetName();
94 G4ThreeVector pdirection =
95 fastTrack.GetPrimaryTrack()->GetMomentum().unit();
96 G4ThreeVector localdir = fastTrack.GetPrimaryTrackLocalDirection();
98 G4ThreeVector worldPosition = fastTrack.GetPrimaryTrack()->GetPosition();
99 G4ThreeVector localPosition = fastTrack.GetPrimaryTrackLocalPosition();
101 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
102 G4double globalTime = fastTrack.GetPrimaryTrack()->GetGlobalTime();
104 G4String particleName =
105 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
107 fastStep.KillPrimaryTrack();
108 fastStep.SetPrimaryTrackPathLength(0.0);
110 if (particleName ==
"kaon+") {
112 }
else if (particleName ==
"kaon-") {
114 }
else if (particleName ==
"anti_proton") {
115 particleName =
"anti-proton";
118 fGarfieldPhysics->
DoIt(particleName, ekin_MeV, globalTime,
119 localPosition.x() / CLHEP::cm, localPosition.y() / CLHEP::cm,
120 localPosition.z() / CLHEP::cm, localdir.x(), localdir.y(),
126 std::vector<GarfieldParticle*>* secondaryParticles =
129 if (secondaryParticles->size() > 0) {
130 fastStep.SetNumberOfSecondaryTracks(secondaryParticles->size());
132 G4double totalEnergySecondaries_MeV = 0;
134 for (std::vector<GarfieldParticle*>::iterator it =
135 secondaryParticles->begin();
136 it != secondaryParticles->end(); ++it) {
137 G4double x = (*it)->getX_mm();
138 G4double y = (*it)->getY_mm();
139 G4double z = (*it)->getZ_mm();
140 G4double eKin_MeV = (*it)->getEkin_MeV();
141 G4double dx = (*it)->getDX();
142 G4double dy = (*it)->getDY();
143 G4double dz = (*it)->getDZ();
144 G4double time = (*it)->getTime();
145 G4ThreeVector momentumDirection(dx, dy, dz);
146 G4ThreeVector position(x, y, z);
147 if ((*it)->getParticleName() ==
"e-") {
148 G4DynamicParticle particle(G4Electron::ElectronDefinition(),
149 momentumDirection, eKin_MeV);
150 fastStep.CreateSecondaryTrack(particle, position, time,
152 totalEnergySecondaries_MeV += eKin_MeV;
153 }
else if ((*it)->getParticleName() ==
"gamma") {
154 G4DynamicParticle particle(G4Gamma::GammaDefinition(),
155 momentumDirection, eKin_MeV);
156 fastStep.CreateSecondaryTrack(particle, position, time,
158 totalEnergySecondaries_MeV += eKin_MeV;
virtual G4bool ModelTrigger(const G4FastTrack &)
~GarfieldG4FastSimulationModel()
virtual void DoIt(const G4FastTrack &, G4FastStep &)
void WriteGeometryToGDML(G4VPhysicalVolume *physicalVolume)
GarfieldG4FastSimulationModel(G4String, G4Region *)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
bool GetCreateSecondariesInGeant4()
static GarfieldPhysics * GetInstance()
std::vector< GarfieldParticle * > * GetSecondaryParticles()
double GetEnergyDeposit_MeV()
void DoIt(std::string particleName, double ekin_MeV, double time, double x_cm, double y_cm, double z_cm, double dx, double dy, double dz)
bool FindParticleNameEnergy(std::string name, double ekin_MeV, std::string program="garfield")
bool FindParticleName(const std::string name, std::string program="garfield")