Garfield++ 4.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 41 of file GarfieldG4FastSimulationModel.cc.

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

◆ GarfieldG4FastSimulationModel() [2/2]

GarfieldG4FastSimulationModel::GarfieldG4FastSimulationModel ( G4String  modelName)

Definition at line 48 of file GarfieldG4FastSimulationModel.cc.

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

◆ ~GarfieldG4FastSimulationModel()

GarfieldG4FastSimulationModel::~GarfieldG4FastSimulationModel ( )

Definition at line 54 of file GarfieldG4FastSimulationModel.cc.

54{}

Member Function Documentation

◆ DoIt()

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

Definition at line 85 of file GarfieldG4FastSimulationModel.cc.

86 {
87
88 G4ThreeVector localdir = fastTrack.GetPrimaryTrackLocalDirection();
89 G4ThreeVector localpos = fastTrack.GetPrimaryTrackLocalPosition();
90
91 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
92 double globalTime = fastTrack.GetPrimaryTrack()->GetGlobalTime();
93
94 G4String particleName =
95 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
96
97 fastStep.KillPrimaryTrack();
98 fastStep.SetPrimaryTrackPathLength(0.0);
99
100 if (particleName == "kaon+") {
101 particleName = "K+";
102 } else if (particleName == "kaon-") {
103 particleName = "K-";
104 } else if (particleName == "anti_proton") {
105 particleName = "anti-proton";
106 }
107
108 fGarfieldPhysics->DoIt(
109 particleName, ekin_MeV, globalTime, localpos.x() / CLHEP::cm,
110 localpos.y() / CLHEP::cm, localpos.z() / CLHEP::cm,
111 localdir.x(), localdir.y(), localdir.z());
112
113 fastStep.SetTotalEnergyDeposited(fGarfieldPhysics->GetEnergyDeposit_MeV());
114
115 if (!fGarfieldPhysics->GetCreateSecondariesInGeant4()) return;
116 std::vector<GarfieldParticle*>* secondaryParticles =
117 fGarfieldPhysics->GetSecondaryParticles();
118
119 if (secondaryParticles->empty()) return;
120 fastStep.SetNumberOfSecondaryTracks(secondaryParticles->size());
121
122 G4double totalEnergySecondaries_MeV = 0;
123
124 for (auto it = secondaryParticles->begin(); it != secondaryParticles->end(); ++it) {
125 G4double eKin_MeV = (*it)->getEkin_MeV();
126 G4double time = (*it)->getTime();
127 G4ThreeVector momentumDirection((*it)->getDX(), (*it)->getDY(),
128 (*it)->getDZ());
129 G4ThreeVector position((*it)->getX_mm(), (*it)->getY_mm(),
130 (*it)->getZ_mm());
131 if ((*it)->getParticleName() == "e-") {
132 G4DynamicParticle particle(G4Electron::ElectronDefinition(),
133 momentumDirection, eKin_MeV);
134 fastStep.CreateSecondaryTrack(particle, position, time, true);
135 totalEnergySecondaries_MeV += eKin_MeV;
136 } else if ((*it)->getParticleName() == "gamma") {
137 G4DynamicParticle particle(G4Gamma::GammaDefinition(),
138 momentumDirection, eKin_MeV);
139 fastStep.CreateSecondaryTrack(particle, position, time, true);
140 totalEnergySecondaries_MeV += eKin_MeV;
141 }
142 }
143
144}
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 64 of file GarfieldG4FastSimulationModel.cc.

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

◆ ModelTrigger()

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

Definition at line 73 of file GarfieldG4FastSimulationModel.cc.

74 {
75 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
76 G4String particleName =
77 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
78 if (fGarfieldPhysics->FindParticleNameEnergy(particleName, ekin_MeV,
79 "garfield")) {
80 return true;
81 }
82 return false;
83}
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 56 of file GarfieldG4FastSimulationModel.cc.

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

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