VMC Examples Version 6.6
Loading...
Searching...
No Matches
GarfieldPhysics.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Virtual Monte Carlo examples
3// Copyright (C) 2007 - 2016 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
10/// \file ExGarfield/src/GarfieldPhysics.cxx
11/// \brief Definition of the GarfieldPhysics class
12///
13/// Garfield garfieldpp example adapted to Virtual Monte Carlo.
14/// This class is imported from garfieldpp example with no modification
15/// in the code.
16///
17/// \date 28/10/2015
18/// \author D. Pheiffer, CERN
19
20#include "GarfieldPhysics.h"
21
22#include "TGeoBBox.h"
23#include "TGeoManager.h"
24
25// I.H.
26//#include "TGDMLParse.h"
27//#include "G4SystemOfUnits.hh"
28
30//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
31
39
45
47{
49 fSecondaryParticles = new std::vector<GarfieldParticle*>();
51 fSensor = 0;
52 fAvalanche = 0;
53 fDrift = 0;
55 fTrackHeed = 0;
56 fGeoManager = 0;
57 fGeometryRoot = 0;
59 fTube = 0;
61 fIonizationModel = "Geant4";
62 Clear();
63}
64
66{
70 delete fMediumMagboltz;
71 delete fSensor;
72 delete fAvalanche;
73 delete fDrift;
75 delete fTrackHeed;
76 delete fGeometryRoot;
77 delete fGeometrySimple;
78
79 std::cout << "Deconstructor GarfieldPhysics" << std::endl;
80}
81
83
84void GarfieldPhysics::SetIonizationModel(std::string model, bool useDefaults)
85{
86 if (model != "Geant4" && model != "Heed") {
87
88 std::cout << "Unknown ionization model " << model << std::endl;
89 std::cout << "Using Geant4 as default model!" << std::endl;
90 model = "Geant4";
91 }
92 fIonizationModel = model;
93
94 if (fIonizationModel == "Geant4") {
95 if (useDefaults == true) {
96 // Particle types and energies for which the G4FastSimulationModel with
97 // Garfield++ is valid
98 this->AddParticleName("e-", 1e-6, 1e-3);
99 this->AddParticleName("gamma", 1e-6, 1e+8);
100 }
101 }
102 else if (fIonizationModel == "Heed") {
103 if (useDefaults == true) {
104 // Particle types and energies for which the G4FastSimulationModel with
105 // Garfield++ is valid
106 this->AddParticleName("gamma", 1e-6, 1e+8);
107 this->AddParticleName("e-", 6e-2, 1e+7);
108 this->AddParticleName("e+", 6e-2, 1e+7);
109 this->AddParticleName("mu-", 1e+1, 1e+8);
110 this->AddParticleName("mu+", 1e+1, 1e+8);
111 this->AddParticleName("pi-", 2e+1, 1e+8);
112 this->AddParticleName("pi+", 2e+1, 1e+8);
113 this->AddParticleName("kaon-", 1e+1, 1e+8);
114 this->AddParticleName("kaon+", 1e+1, 1e+8);
115 this->AddParticleName("proton", 9.e+1, 1e+8);
116 this->AddParticleName("anti_proton", 9.e+1, 1e+8);
117 this->AddParticleName("deuteron", 2.e+2, 1e+8);
118 this->AddParticleName("alpha", 4.e+2, 1e+8);
119 }
120 }
121}
122
124 const std::string particleName, double ekin_min_MeV, double ekin_max_MeV)
125{
126 if (ekin_min_MeV >= ekin_max_MeV) {
127 std::cout << "Ekin_min=" << ekin_min_MeV
128 << " keV is larger than Ekin_max=" << ekin_max_MeV << " keV"
129 << std::endl;
130 return;
131 }
132 std::cout << "Garfield model is applicable for G4Particle " << particleName
133 << " between " << ekin_min_MeV << " keV and " << ekin_max_MeV
134 << " keV" << std::endl;
135 fMapParticlesEnergy->insert(
136 std::make_pair(particleName, std::make_pair(ekin_min_MeV, ekin_max_MeV)));
137}
138
140{
141 MapParticlesEnergy::iterator it;
142 it = fMapParticlesEnergy->find(name);
143 if (it != fMapParticlesEnergy->end()) {
144 return true;
145 }
146 return false;
147}
148
149bool GarfieldPhysics::FindParticleNameEnergy(std::string name, double ekin_MeV)
150{
151 MapParticlesEnergy::iterator it;
152 it = fMapParticlesEnergy->find(name);
153 if (it != fMapParticlesEnergy->end()) {
154 EnergyRange_MeV range = it->second;
155 if (range.first <= ekin_MeV && range.second >= ekin_MeV) {
156 return true;
157 }
158 }
159 return false;
160}
161
163{
164
165 fMediumMagboltz = new Garfield::MediumMagboltz();
166
167 fMediumMagboltz->SetComposition("ar", 70., "co2", 30.);
168 fMediumMagboltz->SetTemperature(293.15);
169 fMediumMagboltz->SetPressure(760.);
170 fMediumMagboltz->EnableDebugging();
171 fMediumMagboltz->Initialise(true);
172 fMediumMagboltz->DisableDebugging();
173 // Set the Penning transfer efficiency.
174 const double rPenning = 0.57;
175 const double lambdaPenning = 0.;
176 fMediumMagboltz->EnablePenningTransfer(rPenning, lambdaPenning, "ar");
177 fMediumMagboltz->LoadGasFile("ar_70_co2_30_1000mbar.gas");
178
179 fSensor = new Garfield::Sensor();
180 fDrift = new Garfield::AvalancheMC();
181 fAvalanche = new Garfield::AvalancheMicroscopic();
182 fComponentAnalyticField = new Garfield::ComponentAnalyticField();
183
185
186 fDrift->SetSensor(fSensor);
187 fAvalanche->SetSensor(fSensor);
188
189 fTrackHeed = new Garfield::TrackHeed();
190 fTrackHeed->EnableDebugging();
191 fTrackHeed->SetSensor(fSensor);
192
193 fTrackHeed->EnableDeltaElectronTransport();
194}
195
197{
198 // Wire radius [cm]
199 const double rWire = 25.e-4;
200 // Outer radius of the tube [cm]
201 const double rTube = 1.451;
202 // Half-length of the tube [cm]
203 const double lTube = 10.;
204
205 fGeometrySimple = new Garfield::GeometrySimple();
206 // Make a tube (centered at the origin, inner radius: 0, outer radius: rTube).
207 fTube = new Garfield::SolidTube(0., 0., rWire, rTube, lTube);
208 // Add the solid to the geometry, together with the medium inside.
211
212 // Voltages
213 const double vWire = 1000.;
214 const double vTube = 0.;
215 // Add the wire in the center.
216 fComponentAnalyticField->AddWire(0., 0., 2 * rWire, vWire, "w");
217 // Add the tube.
218 fComponentAnalyticField->AddTube(rTube, vTube, 0, "t");
219
220 fSensor->AddComponent(fComponentAnalyticField);
221}
222
223void GarfieldPhysics::DoIt(std::string particleName, double ekin_MeV,
224 double time, double x_cm, double y_cm, double z_cm, double dx, double dy,
225 double dz)
226{
227
229
230 // Wire radius [cm]
231 const double rWire = 25.e-4;
232 // Outer radius of the tube [cm]
233 const double rTube = 1.45;
234 // Half-length of the tube [cm]
235 const double lTube = 10.;
236
237 double eKin_eV = ekin_MeV * 1e+6;
238
239 double xc = 0., yc = 0., zc = 0., tc = 0.;
240 // Number of electrons produced in a collision
241 int nc = 0;
242 // Energy loss in a collision
243 double ec = 0.;
244 // Dummy variable (not used at present)
245 double extra = 0.;
246 fEnergyDeposit = 0;
247 // fAvalancheSize = 0;
248 // fGain = 0;
249
250 if (fIonizationModel != "Heed" || particleName == "gamma") {
251 if (particleName == "gamma") {
252 fTrackHeed->TransportPhoton(
253 x_cm, y_cm, z_cm, time, eKin_eV, dx, dy, dz, nc);
254 }
255 else {
256 fTrackHeed->TransportDeltaElectron(
257 x_cm, y_cm, z_cm, time, eKin_eV, dx, dy, dz, nc);
258 fEnergyDeposit = eKin_eV;
259 }
260
261 for (int cl = 0; cl < nc; cl++) {
262 double xe, ye, ze, te;
263 double ee, dxe, dye, dze;
264 fTrackHeed->GetElectron(cl, xe, ye, ze, te, ee, dxe, dye, dze);
265 if (ze < lTube && ze > -lTube && sqrt(xe * xe + ye * ye) < rTube) {
266 nsum++;
267 if (particleName == "gamma") {
268 fEnergyDeposit += fTrackHeed->GetW();
269 }
270 // TO DO: add this histo in MCApplication
271 // analysisManager->FillH3(1, ze * 10, xe * 10, ye * 10);
273 double newTime = te;
274 if (newTime < time) {
275 newTime += time;
276 }
277 fSecondaryParticles->push_back(
278 new GarfieldParticle("e-", ee, newTime, xe, ye, ze, dxe, dye, dze));
279 }
280
281 fDrift->DriftElectron(xe, ye, ze, te);
282
283 double xe1, ye1, ze1, te1;
284 double xe2, ye2, ze2, te2;
285
286 int status;
287 fDrift->GetElectronEndpoint(
288 0, xe1, ye1, ze1, te1, xe2, ye2, ze2, te2, status);
289
290 if (0 < xe2 && xe2 < rWire) {
291 xe2 += 2 * rWire;
292 }
293 if (0 > xe2 && xe2 > -rWire) {
294 xe2 += -2 * rWire;
295 }
296 if (0 < ye2 && ye2 < rWire) {
297 ye2 += 2 * rWire;
298 }
299 if (0 > ye2 && ye2 > -rWire) {
300 ye2 += -2 * rWire;
301 }
302
303 double e2 = 0.1;
304 fAvalanche->AvalancheElectron(xe2, ye2, ze2, te2, e2, 0, 0, 0);
305
306 int ne = 0, ni = 0;
307 fAvalanche->GetAvalancheSize(ne, ni);
308 fAvalancheSize += ne;
309 }
310 }
311 }
312 else {
313 fTrackHeed->SetParticle(particleName);
314 fTrackHeed->SetKineticEnergy(eKin_eV);
315 fTrackHeed->NewTrack(x_cm, y_cm, z_cm, time, dx, dy, dz);
316
317 while (fTrackHeed->GetCluster(xc, yc, zc, tc, nc, ec, extra)) {
318 if (zc < lTube && zc > -lTube && sqrt(xc * xc + yc * yc) < rTube) {
319 nsum += nc;
320 fEnergyDeposit += ec;
321 for (int cl = 0; cl < nc; cl++) {
322 double xe, ye, ze, te;
323 double ee, dxe, dye, dze;
324 fTrackHeed->GetElectron(cl, xe, ye, ze, te, ee, dxe, dye, dze);
325 if (ze < lTube && ze > -lTube && sqrt(xe * xe + ye * ye) < rTube) {
326 // TO DO: add this histo in MCApplication
327 // analysisManager->FillH3(1, ze * 10, xe * 10, ye * 10);
329 double newTime = te;
330 if (newTime < time) {
331 newTime += time;
332 }
334 "e-", ee, newTime, xe, ye, ze, dxe, dye, dze));
335 }
336
337 fDrift->DriftElectron(xe, ye, ze, te);
338
339 double xe1, ye1, ze1, te1;
340 double xe2, ye2, ze2, te2;
341
342 int status;
343 fDrift->GetElectronEndpoint(
344 0, xe1, ye1, ze1, te1, xe2, ye2, ze2, te2, status);
345
346 if (0 < xe2 && xe2 < rWire) {
347 xe2 += 2 * rWire;
348 }
349 if (0 > xe2 && xe2 > -rWire) {
350 xe2 += -2 * rWire;
351 }
352 if (0 < ye2 && ye2 < rWire) {
353 ye2 += 2 * rWire;
354 }
355 if (0 > ye2 && ye2 > -rWire) {
356 ye2 += -2 * rWire;
357 }
358
359 double e2 = 0.1;
360 fAvalanche->AvalancheElectron(xe2, ye2, ze2, te2, e2, 0, 0, 0);
361
362 int ne = 0, ni = 0;
363 fAvalanche->GetAvalancheSize(ne, ni);
364 fAvalancheSize += ne;
365 }
366 }
367 }
368 }
369 }
370 if (nsum > 0) {
372 }
373}
374
375std::vector<GarfieldParticle*>* GarfieldPhysics::GetSecondaryParticles()
376{
377 return fSecondaryParticles;
378}
379
381{
382 if (!fSecondaryParticles->empty()) {
383 fSecondaryParticles->erase(
385 }
386}
Definition of the GarfieldPhysics class.
std::pair< double, double > EnergyRange_MeV
std::map< const std::string, EnergyRange_MeV > MapParticlesEnergy
std::vector< GarfieldParticle * > * fSecondaryParticles
static void Dispose()
Garfield::AvalancheMC * fDrift
static GarfieldPhysics * GetInstance()
Garfield::AvalancheMicroscopic * fAvalanche
TGeoManager * fGeoManager
std::vector< GarfieldParticle * > * GetSecondaryParticles()
void SetIonizationModel(std::string model, bool useDefaults=true)
Garfield::Sensor * fSensor
Garfield::GeometrySimple * fGeometrySimple
Garfield::TrackHeed * fTrackHeed
MapParticlesEnergy * fMapParticlesEnergy
Garfield::GeometryRoot * fGeometryRoot
std::string fIonizationModel
std::string GetIonizationModel()
bool FindParticleNameEnergy(std::string name, double ekin_keV)
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)
Garfield::SolidTube * fTube
void AddParticleName(const std::string particleName, double ekin_min_keV, double ekin_max_keV)
static GarfieldPhysics * fGarfieldPhysics
Garfield::MediumMagboltz * fMediumMagboltz
Garfield::ComponentAnalyticField * fComponentAnalyticField
bool FindParticleName(const std::string name)