Garfield++ 5.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
GarfieldG4FastSimulationModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26/// \file GarfieldG4FastSimulationModel.cc
27/// \brief Implementation of the GarfieldG4FastSimulationModel class
28
30
31#include <iostream>
32
33#include "G4Electron.hh"
34#include "G4GDMLParser.hh"
35#include "G4Gamma.hh"
36#include "G4SystemOfUnits.hh"
37#include "G4VPhysicalVolume.hh"
38
40 G4Region* envelope)
41 : G4VFastSimulationModel(modelName, envelope) {
42 fGarfieldPhysics = GarfieldPhysics::GetInstance();
43 fGarfieldPhysics->InitializePhysics();
44}
45
47 : G4VFastSimulationModel(modelName) {
48 fGarfieldPhysics = GarfieldPhysics::GetInstance();
49 fGarfieldPhysics->InitializePhysics();
50}
51
53
55 G4VPhysicalVolume* physicalVolume) {
56 G4GDMLParser* parser = new G4GDMLParser();
57 remove("garfieldGeometry.gdml");
58 parser->Write("garfieldGeometry.gdml", physicalVolume, false);
59 delete parser;
60}
61
63 const G4ParticleDefinition& particleType) {
64 G4String particleName = particleType.GetParticleName();
65 if (fGarfieldPhysics->FindParticleName(particleName, "garfield")) {
66 return true;
67 }
68 return false;
69}
70
72 const G4FastTrack& fastTrack) {
73 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
74 G4String particleName =
75 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
76 if (fGarfieldPhysics->FindParticleNameEnergy(particleName, ekin_MeV,
77 "garfield")) {
78 return true;
79 }
80 return false;
81}
82
83void GarfieldG4FastSimulationModel::DoIt(const G4FastTrack& fastTrack,
84 G4FastStep& fastStep) {
85
86 G4ThreeVector localdir = fastTrack.GetPrimaryTrackLocalDirection();
87 G4ThreeVector localpos = fastTrack.GetPrimaryTrackLocalPosition();
88
89 double ekin_MeV = fastTrack.GetPrimaryTrack()->GetKineticEnergy() / MeV;
90 double globalTime = fastTrack.GetPrimaryTrack()->GetGlobalTime();
91
92 G4String particleName =
93 fastTrack.GetPrimaryTrack()->GetParticleDefinition()->GetParticleName();
94
95 fastStep.KillPrimaryTrack();
96 fastStep.SetPrimaryTrackPathLength(0.0);
97
98 if (particleName == "kaon+") {
99 particleName = "K+";
100 } else if (particleName == "kaon-") {
101 particleName = "K-";
102 } else if (particleName == "anti_proton") {
103 particleName = "anti-proton";
104 }
105
106 fGarfieldPhysics->DoIt(
107 particleName, ekin_MeV, globalTime, localpos.x() / CLHEP::cm,
108 localpos.y() / CLHEP::cm, localpos.z() / CLHEP::cm,
109 localdir.x(), localdir.y(), localdir.z());
110
111 fastStep.SetTotalEnergyDeposited(fGarfieldPhysics->GetEnergyDeposit_MeV());
112
113 if (!fGarfieldPhysics->GetCreateSecondariesInGeant4()) return;
114 const auto& secondaryParticles = fGarfieldPhysics->GetSecondaryParticles();
115
116 if (secondaryParticles.empty()) return;
117 fastStep.SetNumberOfSecondaryTracks(secondaryParticles.size());
118
119 G4double totalEnergySecondaries_MeV = 0.;
120
121 for (const auto& sp : secondaryParticles) {
122 G4double eKin_MeV = sp.getEkin_MeV();
123 G4double time = sp.getTime();
124 G4ThreeVector momentumDirection(sp.getDX(), sp.getDY(), sp.getDZ());
125 G4ThreeVector position(sp.getX_mm(), sp.getY_mm(), sp.getZ_mm());
126 if (sp.getParticleName() == "e-") {
127 G4DynamicParticle particle(G4Electron::ElectronDefinition(),
128 momentumDirection, eKin_MeV);
129 fastStep.CreateSecondaryTrack(particle, position, time, true);
130 totalEnergySecondaries_MeV += eKin_MeV;
131 } else if (sp.getParticleName() == "gamma") {
132 G4DynamicParticle particle(G4Gamma::GammaDefinition(),
133 momentumDirection, eKin_MeV);
134 fastStep.CreateSecondaryTrack(particle, position, time, true);
135 totalEnergySecondaries_MeV += eKin_MeV;
136 }
137 }
138
139}
virtual G4bool ModelTrigger(const G4FastTrack &)
virtual void DoIt(const G4FastTrack &, G4FastStep &)
void WriteGeometryToGDML(G4VPhysicalVolume *physicalVolume)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
static GarfieldPhysics * GetInstance()