VMC Examples Version 6.6
Loading...
Searching...
No Matches
GarfieldG4FastSimulationModel Class Reference

#include <GarfieldG4FastSimulationModel.h>

Inheritance diagram for GarfieldG4FastSimulationModel:

Public Member Functions

 GarfieldG4FastSimulationModel (G4String, G4Region *)
 
 GarfieldG4FastSimulationModel (G4String)
 
 ~GarfieldG4FastSimulationModel ()
 
void SetPhysics (GarfieldPhysics *fGarfieldPhysics)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual G4bool ModelTrigger (const G4FastTrack &)
 
virtual void DoIt (const G4FastTrack &, G4FastStep &)
 

Private Attributes

GarfieldPhysicsfGarfieldPhysics
 

Detailed Description

Definition at line 36 of file GarfieldG4FastSimulationModel.h.

Constructor & Destructor Documentation

◆ GarfieldG4FastSimulationModel() [1/2]

GarfieldG4FastSimulationModel::GarfieldG4FastSimulationModel ( G4String modelName,
G4Region * envelope )

◆ GarfieldG4FastSimulationModel() [2/2]

GarfieldG4FastSimulationModel::GarfieldG4FastSimulationModel ( G4String modelName)

◆ ~GarfieldG4FastSimulationModel()

GarfieldG4FastSimulationModel::~GarfieldG4FastSimulationModel ( )

Definition at line 54 of file GarfieldG4FastSimulationModel.cxx.

54{}

Member Function Documentation

◆ SetPhysics()

void GarfieldG4FastSimulationModel::SetPhysics ( GarfieldPhysics * fGarfieldPhysics)

◆ IsApplicable()

G4bool GarfieldG4FastSimulationModel::IsApplicable ( const G4ParticleDefinition & particleType)
virtual

Definition at line 66 of file GarfieldG4FastSimulationModel.cxx.

68{
69
70 G4String particleName = particleType.GetParticleName();
71
72 if (fGarfieldPhysics->FindParticleName(particleName)) {
73 return true;
74 }
75 return false;
76}
bool FindParticleName(const std::string name)

◆ ModelTrigger()

G4bool GarfieldG4FastSimulationModel::ModelTrigger ( const G4FastTrack & fastTrack)
virtual

Definition at line 78 of file GarfieldG4FastSimulationModel.cxx.

79{
80 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
81 G4String particleName =
82 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
83 if (fGarfieldPhysics->FindParticleNameEnergy(particleName, ekin_MeV)) {
84 return true;
85 }
86 return false;
87}
bool FindParticleNameEnergy(std::string name, double ekin_keV)

◆ DoIt()

void GarfieldG4FastSimulationModel::DoIt ( const G4FastTrack & fastTrack,
G4FastStep & fastStep )
virtual

Definition at line 89 of file GarfieldG4FastSimulationModel.cxx.

91{
92
93 G4TouchableHandle theTouchable =
94 fastTrack.GetPrimaryTrack()->GetTouchableHandle();
95 G4String name = theTouchable->GetVolume()->GetName();
96
97 G4ThreeVector pdirection = fastTrack.GetPrimaryTrack()->GetMomentum().unit();
98 G4ThreeVector localdir = fastTrack.GetPrimaryTrackLocalDirection();
99
100 G4ThreeVector worldPosition = fastTrack.GetPrimaryTrack()->GetPosition();
101 G4ThreeVector localPosition = fastTrack.GetPrimaryTrackLocalPosition();
102
103 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
104 G4double globalTime = fastTrack.GetPrimaryTrack()->GetGlobalTime();
105
106 G4String particleName =
107 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
108
109 fastStep.KillPrimaryTrack();
110 fastStep.SetPrimaryTrackPathLength(0.0);
111
112 if (particleName == "kaon+") {
113 particleName = "K+";
114 }
115 else if (particleName == "kaon-") {
116 particleName = "K-";
117 }
118 else if (particleName == "anti_proton") {
119 particleName = "anti-proton";
120 }
121
122 fGarfieldPhysics->DoIt(particleName, ekin_MeV, globalTime,
123 localPosition.x() / CLHEP::cm, localPosition.y() / CLHEP::cm,
124 localPosition.z() / CLHEP::cm, localdir.x(), localdir.y(), localdir.z());
125
126 fastStep.SetTotalEnergyDeposited(fGarfieldPhysics->GetEnergyDeposit_MeV());
127
129 std::vector<GarfieldParticle*>* secondaryParticles =
131
132 if (secondaryParticles->size() > 0) {
133 fastStep.SetNumberOfSecondaryTracks(secondaryParticles->size());
134
135 G4double totalEnergySecondaries_MeV = 0;
136
137 for (std::vector<GarfieldParticle*>::iterator it =
138 secondaryParticles->begin();
139 it != secondaryParticles->end(); ++it) {
140 G4double x = (*it)->getX_mm();
141 G4double y = (*it)->getY_mm();
142 G4double z = (*it)->getZ_mm();
143 G4double eKin_MeV = (*it)->getEkin_MeV();
144 G4double dx = (*it)->getDX();
145 G4double dy = (*it)->getDY();
146 G4double dz = (*it)->getDZ();
147 G4double time = (*it)->getTime();
148 G4ThreeVector momentumDirection(dx, dy, dz);
149 G4ThreeVector position(x, y, z);
150 if ((*it)->getParticleName() == "e-") {
151 G4DynamicParticle particle(
152 G4Electron::ElectronDefinition(), momentumDirection, eKin_MeV);
153 fastStep.CreateSecondaryTrack(particle, position, time, true);
154 totalEnergySecondaries_MeV += eKin_MeV;
155 }
156 else if ((*it)->getParticleName() == "gamma") {
157 G4DynamicParticle particle(
158 G4Gamma::GammaDefinition(), momentumDirection, eKin_MeV);
159 fastStep.CreateSecondaryTrack(particle, position, time, true);
160 totalEnergySecondaries_MeV += eKin_MeV;
161 }
162 }
163 }
164 }
165}
bool GetCreateSecondariesInGeant4()
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)

Member Data Documentation

◆ fGarfieldPhysics

GarfieldPhysics* GarfieldG4FastSimulationModel::fGarfieldPhysics
private

Definition at line 55 of file GarfieldG4FastSimulationModel.h.


The documentation for this class was generated from the following files: