Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
GarfieldG4FastSimulationModel Class Reference

#include <GarfieldG4FastSimulationModel.hh>

+ Inheritance diagram for GarfieldG4FastSimulationModel:

Public Member Functions

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

Detailed Description

Definition at line 39 of file GarfieldG4FastSimulationModel.hh.

Constructor & Destructor Documentation

◆ GarfieldG4FastSimulationModel() [1/2]

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

Definition at line 40 of file GarfieldG4FastSimulationModel.cc.

41 :
42 G4VFastSimulationModel(modelName, envelope) {
43 fGarfieldPhysics = GarfieldPhysics::GetInstance();
44 fGarfieldPhysics->InitializePhysics();
45
46}
static GarfieldPhysics * GetInstance()

◆ GarfieldG4FastSimulationModel() [2/2]

GarfieldG4FastSimulationModel::GarfieldG4FastSimulationModel ( G4String  modelName)

Definition at line 48 of file GarfieldG4FastSimulationModel.cc.

48 :
49 G4VFastSimulationModel(modelName) {
50 fGarfieldPhysics = GarfieldPhysics::GetInstance();
51 fGarfieldPhysics->InitializePhysics();
52}

◆ ~GarfieldG4FastSimulationModel()

GarfieldG4FastSimulationModel::~GarfieldG4FastSimulationModel ( )

Definition at line 54 of file GarfieldG4FastSimulationModel.cc.

54 {
55}

Member Function Documentation

◆ DoIt()

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

Definition at line 87 of file GarfieldG4FastSimulationModel.cc.

88 {
89
90 G4TouchableHandle theTouchable =
91 fastTrack.GetPrimaryTrack()->GetTouchableHandle();
92 G4String name = theTouchable->GetVolume()->GetName();
93
94 G4ThreeVector pdirection =
95 fastTrack.GetPrimaryTrack()->GetMomentum().unit();
96 G4ThreeVector localdir = fastTrack.GetPrimaryTrackLocalDirection();
97
98 G4ThreeVector worldPosition = fastTrack.GetPrimaryTrack()->GetPosition();
99 G4ThreeVector localPosition = fastTrack.GetPrimaryTrackLocalPosition();
100
101 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
102 G4double globalTime = fastTrack.GetPrimaryTrack()->GetGlobalTime();
103
104 G4String particleName =
105 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
106
107 fastStep.KillPrimaryTrack();
108 fastStep.SetPrimaryTrackPathLength(0.0);
109
110 if (particleName == "kaon+") {
111 particleName = "K+";
112 } else if (particleName == "kaon-") {
113 particleName = "K-";
114 } else if (particleName == "anti_proton") {
115 particleName = "anti-proton";
116 }
117
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(),
121 localdir.z());
122
123 fastStep.SetTotalEnergyDeposited(fGarfieldPhysics->GetEnergyDeposit_MeV());
124
125 if (fGarfieldPhysics->GetCreateSecondariesInGeant4()) {
126 std::vector<GarfieldParticle*>* secondaryParticles =
127 fGarfieldPhysics->GetSecondaryParticles();
128
129 if (secondaryParticles->size() > 0) {
130 fastStep.SetNumberOfSecondaryTracks(secondaryParticles->size());
131
132 G4double totalEnergySecondaries_MeV = 0;
133
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,
151 true);
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,
157 true);
158 totalEnergySecondaries_MeV += eKin_MeV;
159 }
160
161 }
162 }
163 }
164}
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)

◆ IsApplicable()

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

Definition at line 66 of file GarfieldG4FastSimulationModel.cc.

67 {
68 G4String particleName = particleType.GetParticleName();
69 if (fGarfieldPhysics->FindParticleName(particleName, "garfield")) {
70 return true;
71 }
72 return false;
73}
bool FindParticleName(const std::string name, std::string program="garfield")

◆ ModelTrigger()

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

Definition at line 75 of file GarfieldG4FastSimulationModel.cc.

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

◆ SetPhysics()

void GarfieldG4FastSimulationModel::SetPhysics ( GarfieldPhysics fGarfieldPhysics)

◆ WriteGeometryToGDML()

void GarfieldG4FastSimulationModel::WriteGeometryToGDML ( G4VPhysicalVolume *  physicalVolume)

Definition at line 57 of file GarfieldG4FastSimulationModel.cc.

58 {
59
60 G4GDMLParser* parser = new G4GDMLParser();
61 remove("garfieldGeometry.gdml");
62 parser->Write("garfieldGeometry.gdml", physicalVolume, false);
63 delete parser;
64}

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