Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
GarfieldPhysics.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 GarfieldPhysics.cc
27/// \brief Implementation of the GarfieldPhysics class
28#include "GarfieldPhysics.hh"
29#include "GarfieldAnalysis.hh"
30
33
34GarfieldPhysics* GarfieldPhysics::fGarfieldPhysics = 0;
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
36
38 if (!fGarfieldPhysics) {
39 fGarfieldPhysics = new GarfieldPhysics();
40 }
41 return fGarfieldPhysics;
42}
43
45 delete fGarfieldPhysics;
46 fGarfieldPhysics = 0;
47}
48
49GarfieldPhysics::GarfieldPhysics() {
50 fSecondaryParticles = new std::vector<GarfieldParticle*>();
51 fMediumMagboltz = 0;
52 fSensor = 0;
53 fComponentAnalyticField = 0;
54 fTrackHeed = 0;
55 createSecondariesInGeant4 = false;
56 fIonizationModel = "PAIPhot";
57}
58
59GarfieldPhysics::~GarfieldPhysics() {
61 delete fSecondaryParticles;
62 delete fMediumMagboltz;
63 delete fSensor;
64 delete fComponentAnalyticField;
65 delete fTrackHeed;
66
67 std::cout << "Deconstructor GarfieldPhysics" << std::endl;
68}
69
70std::string GarfieldPhysics::GetIonizationModel() { return fIonizationModel; }
71
72void GarfieldPhysics::SetIonizationModel(std::string model, bool useDefaults) {
73 if (model != "PAIPhot" && model != "PAI" && model != "Heed") {
74 std::cout << "Unknown ionization model " << model << std::endl;
75 std::cout << "Using PAIPhot as default model!" << std::endl;
76 model = "PAIPhot";
77 }
78 fIonizationModel = model;
79
80 if (fIonizationModel == "PAIPhot" || fIonizationModel == "PAI") {
81 if (useDefaults == true) {
82 // Particle types and energies for which the G4FastSimulationModel with
83 // Garfield++ is valid
84 this->AddParticleName("e-", 1e-6, 1e-3, "garfield");
85 this->AddParticleName("gamma", 1e-6, 1e+8, "garfield");
86
87 // Particle types and energies for which the PAI or PAIPhot model is valid
88 this->AddParticleName("e-", 0, 1e+8, "geant4");
89 this->AddParticleName("e+", 0, 1e+8, "geant4");
90 this->AddParticleName("mu-", 0, 1e+8, "geant4");
91 this->AddParticleName("mu+", 0, 1e+8, "geant4");
92 this->AddParticleName("proton", 0, 1e+8, "geant4");
93 this->AddParticleName("pi+", 0, 1e+8, "geant4");
94 this->AddParticleName("pi-", 0, 1e+8, "geant4");
95 this->AddParticleName("alpha", 0, 1e+8, "geant4");
96 this->AddParticleName("He3", 0, 1e+8, "geant4");
97 this->AddParticleName("GenericIon", 0, 1e+8, "geant4");
98 }
99
100 } else if (fIonizationModel == "Heed") {
101 if (useDefaults == true) {
102 // Particle types and energies for which the G4FastSimulationModel with
103 // Garfield++ is valid
104 this->AddParticleName("gamma", 1e-6, 1e+8, "garfield");
105 this->AddParticleName("e-", 6e-2, 1e+7, "garfield");
106 this->AddParticleName("e+", 6e-2, 1e+7, "garfield");
107 this->AddParticleName("mu-", 1e+1, 1e+8, "garfield");
108 this->AddParticleName("mu+", 1e+1, 1e+8, "garfield");
109 this->AddParticleName("pi-", 2e+1, 1e+8, "garfield");
110 this->AddParticleName("pi+", 2e+1, 1e+8, "garfield");
111 this->AddParticleName("kaon-", 1e+1, 1e+8, "garfield");
112 this->AddParticleName("kaon+", 1e+1, 1e+8, "garfield");
113 this->AddParticleName("proton", 9.e+1, 1e+8, "garfield");
114 this->AddParticleName("anti_proton", 9.e+1, 1e+8, "garfield");
115 this->AddParticleName("deuteron", 2.e+2, 1e+8, "garfield");
116 this->AddParticleName("alpha", 4.e+2, 1e+8, "garfield");
117 }
118 }
119}
120
121void GarfieldPhysics::AddParticleName(const std::string particleName,
122 double ekin_min_MeV, double ekin_max_MeV,
123 std::string program) {
124 if (ekin_min_MeV >= ekin_max_MeV) {
125 std::cout << "Ekin_min=" << ekin_min_MeV
126 << " keV is larger than Ekin_max=" << ekin_max_MeV << " keV"
127 << std::endl;
128 return;
129 }
130
131 if (program == "garfield") {
132 std::cout << "Garfield model (Heed) is applicable for G4Particle "
133 << particleName << " between " << ekin_min_MeV << " MeV and "
134 << ekin_max_MeV << " MeV" << std::endl;
135
136 fMapParticlesEnergyGarfield.insert(std::make_pair(
137 particleName, std::make_pair(ekin_min_MeV, ekin_max_MeV)));
138 } else {
139 std::cout << fIonizationModel << " is applicable for G4Particle "
140 << particleName << " between " << ekin_min_MeV << " MeV and "
141 << ekin_max_MeV << " MeV" << std::endl;
142 fMapParticlesEnergyGeant4.insert(std::make_pair(
143 particleName, std::make_pair(ekin_min_MeV, ekin_max_MeV)));
144 }
145}
146
147bool GarfieldPhysics::FindParticleName(std::string name, std::string program) {
148 if (program == "garfield") {
149 auto it = fMapParticlesEnergyGarfield.find(name);
150 if (it != fMapParticlesEnergyGarfield.end()) return true;
151 } else {
152 auto it = fMapParticlesEnergyGeant4.find(name);
153 if (it != fMapParticlesEnergyGeant4.end()) return true;
154 }
155 return false;
156}
157
158bool GarfieldPhysics::FindParticleNameEnergy(std::string name, double ekin_MeV,
159 std::string program) {
160 if (program == "garfield") {
161 auto it = fMapParticlesEnergyGarfield.find(name);
162 if (it != fMapParticlesEnergyGarfield.end()) {
163 EnergyRange_MeV range = it->second;
164 if (range.first <= ekin_MeV && range.second >= ekin_MeV) {
165 return true;
166 }
167 }
168 } else {
169 auto it = fMapParticlesEnergyGeant4.find(name);
170 if (it != fMapParticlesEnergyGeant4.end()) {
171 EnergyRange_MeV range = it->second;
172 if (range.first <= ekin_MeV && range.second >= ekin_MeV) {
173 return true;
174 }
175 }
176 }
177 return false;
178}
179
181 std::string program) {
182 if (program == "garfield") {
183 auto it = fMapParticlesEnergyGarfield.find(name);
184 if (it != fMapParticlesEnergyGarfield.end()) {
185 EnergyRange_MeV range = it->second;
186 return range.first;
187 }
188 } else {
189 auto it = fMapParticlesEnergyGeant4.find(name);
190 if (it != fMapParticlesEnergyGeant4.end()) {
191 EnergyRange_MeV range = it->second;
192 return range.first;
193 }
194 }
195 return -1;
196}
197
199 std::string program) {
200 if (program == "garfield") {
201 auto it = fMapParticlesEnergyGarfield.find(name);
202 if (it != fMapParticlesEnergyGarfield.end()) {
203 EnergyRange_MeV range = it->second;
204 return range.second;
205 }
206 } else {
207 auto it = fMapParticlesEnergyGeant4.find(name);
208 if (it != fMapParticlesEnergyGeant4.end()) {
209 EnergyRange_MeV range = it->second;
210 return range.second;
211 }
212 }
213 return -1;
214}
215
217 // Define the gas mixture.
218 fMediumMagboltz = new Garfield::MediumMagboltz();
219 fMediumMagboltz->SetComposition("ar", 70., "co2", 30.);
220 fMediumMagboltz->SetTemperature(293.15);
221 fMediumMagboltz->SetPressure(760.);
222 fMediumMagboltz->Initialise(true);
223 // Set the Penning transfer efficiency.
224 const double rPenning = 0.57;
225 const double lambdaPenning = 0.;
226 fMediumMagboltz->EnablePenningTransfer(rPenning, lambdaPenning, "ar");
227 fMediumMagboltz->LoadGasFile("ar_70_co2_30_1000mbar.gas");
228
229 fComponentAnalyticField = new Garfield::ComponentAnalyticField();
230 fComponentAnalyticField->SetMedium(fMediumMagboltz);
231 // Wire radius [cm]
232 const double rWire = 25.e-4;
233 // Tube radius [cm]
234 const double rTube = 1.451;
235 // Voltages
236 const double vWire = 1000.;
237 const double vTube = 0.;
238 // Add the wire in the center.
239 fComponentAnalyticField->AddWire(0., 0., 2 * rWire, vWire, "w");
240 // Add the tube.
241 fComponentAnalyticField->AddTube(rTube, vTube, 0, "t");
242
243 fSensor = new Garfield::Sensor();
244 fSensor->AddComponent(fComponentAnalyticField);
245
246 fTrackHeed = new Garfield::TrackHeed();
247 fTrackHeed->SetSensor(fSensor);
248 fTrackHeed->EnableDeltaElectronTransport();
249}
250
251void GarfieldPhysics::DoIt(std::string particleName, double ekin_MeV,
252 double time, double x_cm, double y_cm, double z_cm,
253 double dx, double dy, double dz) {
254 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
255 fEnergyDeposit = 0;
257
259 drift.SetSensor(fSensor);
260 drift.SetDistanceSteps(1.e-4);
262 avalanche.SetSensor(fSensor);
263
264 // Wire radius [cm]
265 const double rWire = 25.e-4;
266 // Outer radius of the tube [cm]
267 const double rTube = 1.45;
268 // Half-length of the tube [cm]
269 const double lTube = 10.;
270
271 double eKin_eV = ekin_MeV * 1e+6;
272
273 double xc = 0., yc = 0., zc = 0., tc = 0.;
274 // Number of electrons produced in a collision
275 int nc = 0;
276 // Energy loss in a collision
277 double ec = 0.;
278 // Dummy variable (not used at present)
279 double extra = 0.;
280 fEnergyDeposit = 0;
281 if (fIonizationModel != "Heed" || particleName == "gamma") {
282 if (particleName == "gamma") {
283 fTrackHeed->TransportPhoton(x_cm, y_cm, z_cm, time, eKin_eV, dx, dy, dz,
284 nc);
285 } else {
286 fTrackHeed->TransportDeltaElectron(x_cm, y_cm, z_cm, time, eKin_eV, dx,
287 dy, dz, nc);
288 fEnergyDeposit = eKin_eV;
289 }
290
291 for (int cl = 0; cl < nc; cl++) {
292 double xe, ye, ze, te;
293 double ee, dxe, dye, dze;
294 fTrackHeed->GetElectron(cl, xe, ye, ze, te, ee, dxe, dye, dze);
295 if (ze < lTube && ze > -lTube && sqrt(xe * xe + ye * ye) < rTube) {
296 nsum++;
297 if (particleName == "gamma") {
298 fEnergyDeposit += fTrackHeed->GetW();
299 }
300 analysisManager->FillH3(1, ze * 10, xe * 10, ye * 10);
301 if (createSecondariesInGeant4) {
302 double newTime = te;
303 if (newTime < time) {
304 newTime += time;
305 }
306 fSecondaryParticles->push_back(new GarfieldParticle(
307 "e-", ee, newTime, xe, ye, ze, dxe, dye, dze));
308 }
309
310 drift.DriftElectron(xe, ye, ze, te);
311
312 double xe1, ye1, ze1, te1;
313 double xe2, ye2, ze2, te2;
314
315 int status;
316 drift.GetElectronEndpoint(0, xe1, ye1, ze1, te1, xe2, ye2, ze2, te2,
317 status);
318
319 if (0 < xe2 && xe2 < rWire) {
320 xe2 += 2 * rWire;
321 } else if (0 > xe2 && xe2 > -rWire) {
322 xe2 += -2 * rWire;
323 }
324 if (0 < ye2 && ye2 < rWire) {
325 ye2 += 2 * rWire;
326 } else if (0 > ye2 && ye2 > -rWire) {
327 ye2 += -2 * rWire;
328 }
329
330 double e2 = 0.1;
331 avalanche.AvalancheElectron(xe2, ye2, ze2, te2, e2, 0, 0, 0);
332
333 int ne = 0, ni = 0;
334 avalanche.GetAvalancheSize(ne, ni);
335 fAvalancheSize += ne;
336 }
337 }
338 } else {
339 fTrackHeed->SetParticle(particleName);
340 fTrackHeed->SetKineticEnergy(eKin_eV);
341 fTrackHeed->NewTrack(x_cm, y_cm, z_cm, time, dx, dy, dz);
342
343 while (fTrackHeed->GetCluster(xc, yc, zc, tc, nc, ec, extra)) {
344 if (zc < lTube && zc > -lTube && sqrt(xc * xc + yc * yc) < rTube) {
345 nsum += nc;
346 fEnergyDeposit += ec;
347 for (int cl = 0; cl < nc; cl++) {
348 double xe, ye, ze, te;
349 double ee, dxe, dye, dze;
350 fTrackHeed->GetElectron(cl, xe, ye, ze, te, ee, dxe, dye, dze);
351 if (ze < lTube && ze > -lTube && sqrt(xe * xe + ye * ye) < rTube) {
352 analysisManager->FillH3(1, ze * 10, xe * 10, ye * 10);
353 if (createSecondariesInGeant4) {
354 double newTime = te;
355 if (newTime < time) {
356 newTime += time;
357 }
358 fSecondaryParticles->push_back(new GarfieldParticle(
359 "e-", ee, newTime, xe, ye, ze, dxe, dye, dze));
360 }
361
362 drift.DriftElectron(xe, ye, ze, te);
363
364 double xe1, ye1, ze1, te1;
365 double xe2, ye2, ze2, te2;
366
367 int status;
368 drift.GetElectronEndpoint(0, xe1, ye1, ze1, te1, xe2, ye2, ze2,
369 te2, status);
370
371 if (0 < xe2 && xe2 < rWire) {
372 xe2 += 2 * rWire;
373 } else if (0 > xe2 && xe2 > -rWire) {
374 xe2 += -2 * rWire;
375 }
376 if (0 < ye2 && ye2 < rWire) {
377 ye2 += 2 * rWire;
378 } else if (0 > ye2 && ye2 > -rWire) {
379 ye2 += -2 * rWire;
380 }
381
382 double e2 = 0.1;
383 avalanche.AvalancheElectron(xe2, ye2, ze2, te2, e2, 0, 0, 0);
384
385 int ne = 0, ni = 0;
386 avalanche.GetAvalancheSize(ne, ni);
387 fAvalancheSize += ne;
388 }
389 }
390 }
391 }
392 }
393 fGain = fAvalancheSize / nsum;
394}
395
396std::vector<GarfieldParticle*>* GarfieldPhysics::GetSecondaryParticles() {
397 return fSecondaryParticles;
398}
399
401 if (!fSecondaryParticles->empty()) {
402 fSecondaryParticles->erase(fSecondaryParticles->begin(),
403 fSecondaryParticles->end());
404 }
405}
Selection of the analysis technology.
Definition of the GarfieldPhysics class.
std::pair< double, double > EnergyRange_MeV
static void Dispose()
double GetMaxEnergyMeVParticle(std::string name, std::string program="garfield")
void DeleteSecondaryParticles()
static GarfieldPhysics * GetInstance()
std::vector< GarfieldParticle * > * GetSecondaryParticles()
void SetIonizationModel(std::string model, bool useDefaults=true)
double GetMinEnergyMeVParticle(std::string name, std::string program="garfield")
void AddParticleName(const std::string particleName, double ekin_min_MeV, double ekin_max_MeV, std::string program)
std::string GetIonizationModel()
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")
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Simulate the drift line of an electron from a given starting point.
Definition: AvalancheMC.cc:186
void GetElectronEndpoint(const unsigned int i, double &x0, double &y0, double &z0, double &t0, double &x1, double &y1, double &z1, double &t1, int &status) const
Definition: AvalancheMC.cc:166
void SetSensor(Sensor *s)
Set the sensor.
Definition: AvalancheMC.cc:29
void SetDistanceSteps(const double d=0.001)
Use fixed distance steps (default 10 um).
Definition: AvalancheMC.cc:62
Calculate electron drift lines and avalanches using microscopic tracking.
void GetAvalancheSize(int &ne, int &ni) const
Return the number of electrons and ions in the avalanche.
void SetSensor(Sensor *sensor)
Set the sensor.
bool AvalancheElectron(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0=0., const double dy0=0., const double dz0=0.)
Calculate an avalanche initiated by a given electron.
void AddTube(const double radius, const double voltage, const int nEdges, const std::string &label)
Add a tube.
void SetMedium(Medium *medium)
Set the medium inside the cell.
void AddWire(const double x, const double y, const double diameter, const double voltage, const std::string &label, const double length=100., const double tension=50., const double rho=19.3, const int ntrap=5)
Add a wire at (x, y) .
bool SetComposition(const std::string &gas1, const double f1=1., const std::string &gas2="", const double f2=0., const std::string &gas3="", const double f3=0., const std::string &gas4="", const double f4=0., const std::string &gas5="", const double f5=0., const std::string &gas6="", const double f6=0.)
Set the gas mixture.
Definition: MediumGas.cc:135
bool LoadGasFile(const std::string &filename)
Read table of gas properties (transport parameters) from file.
Definition: MediumGas.cc:312
bool EnablePenningTransfer(const double r, const double lambda) override
bool Initialise(const bool verbose=false)
void SetTemperature(const double t)
Set the temperature [K].
Definition: Medium.cc:71
void SetPressure(const double p)
Definition: Medium.cc:81
void AddComponent(Component *comp)
Add a component.
Definition: Sensor.cc:316
Generate tracks using Heed++.
Definition: TrackHeed.hh:35
bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra) override
Definition: TrackHeed.cc:222
void TransportDeltaElectron(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0, const double dy0, const double dz0, int &ne, int &ni)
Definition: TrackHeed.cc:422
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
Definition: TrackHeed.cc:59
void EnableDeltaElectronTransport()
Switch simulation of delta electrons on.
Definition: TrackHeed.hh:159
bool GetElectron(const unsigned int i, double &x, double &y, double &z, double &t, double &e, double &dx, double &dy, double &dz)
Definition: TrackHeed.cc:354
void TransportPhoton(const double x0, const double y0, const double z0, const double t0, const double e0, const double dx0, const double dy0, const double dz0, int &ne, int &ni)
Definition: TrackHeed.cc:535
double GetW() const
Return the W value of the medium (of the last simulated track).
Definition: TrackHeed.cc:1047
void SetSensor(Sensor *s)
Set the sensor through which to transport the particle.
Definition: Track.cc:166
void SetKineticEnergy(const double ekin)
Set the kinetic energy of the particle.
Definition: Track.cc:153
virtual void SetParticle(const std::string &part)
Definition: Track.cc:15