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

#include <GarfieldPhysics.hh>

Public Member Functions

void InitializePhysics ()
 
void CreateGeometry ()
 
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)
 
void AddParticleName (const std::string particleName, double ekin_min_MeV, double ekin_max_MeV, std::string program)
 
bool FindParticleName (const std::string name, std::string program="garfield")
 
bool FindParticleNameEnergy (std::string name, double ekin_MeV, std::string program="garfield")
 
double GetMinEnergyMeVParticle (std::string name, std::string program="garfield")
 
double GetMaxEnergyMeVParticle (std::string name, std::string program="garfield")
 
void SetIonizationModel (std::string model, bool useDefaults=true)
 
std::string GetIonizationModel ()
 
std::vector< GarfieldParticle * > * GetSecondaryParticles ()
 
void DeleteSecondaryParticles ()
 
void EnableCreateSecondariesInGeant4 (bool flag)
 
bool GetCreateSecondariesInGeant4 ()
 
double GetEnergyDeposit_MeV ()
 
double GetAvalancheSize ()
 
double GetGain ()
 
void Clear ()
 

Static Public Member Functions

static GarfieldPhysicsGetInstance ()
 
static void Dispose ()
 

Detailed Description

Definition at line 76 of file GarfieldPhysics.hh.

Member Function Documentation

◆ AddParticleName()

void GarfieldPhysics::AddParticleName ( const std::string  particleName,
double  ekin_min_MeV,
double  ekin_max_MeV,
std::string  program 
)

Definition at line 140 of file GarfieldPhysics.cc.

141 {
142 if (ekin_min_MeV >= ekin_max_MeV) {
143 std::cout << "Ekin_min=" << ekin_min_MeV
144 << " keV is larger than Ekin_max=" << ekin_max_MeV << " keV"
145 << std::endl;
146 return;
147 }
148
149 if (program == "garfield") {
150 std::cout << "Garfield model (Heed) is applicable for G4Particle "
151 << particleName << " between " << ekin_min_MeV << " MeV and "
152 << ekin_max_MeV << " MeV" << std::endl;
153
154 fMapParticlesEnergyGarfield->insert(
155 std::make_pair(particleName,
156 std::make_pair(ekin_min_MeV, ekin_max_MeV)));
157 } else {
158 std::cout << fIonizationModel << " is applicable for G4Particle "
159 << particleName << " between " << ekin_min_MeV << " MeV and "
160 << ekin_max_MeV << " MeV" << std::endl;
161 fMapParticlesEnergyGeant4->insert(
162 std::make_pair(particleName,
163 std::make_pair(ekin_min_MeV, ekin_max_MeV)));
164 }
165
166}

Referenced by SetIonizationModel(), and GarfieldMessenger::SetNewValue().

◆ Clear()

void GarfieldPhysics::Clear ( )
inline

Definition at line 101 of file GarfieldPhysics.hh.

101{fEnergyDeposit=0;fAvalancheSize=0;fGain=0;nsum=0;}

Referenced by GarfieldEventAction::BeginOfEventAction().

◆ CreateGeometry()

void GarfieldPhysics::CreateGeometry ( )

Definition at line 282 of file GarfieldPhysics.cc.

282 {
283// Wire radius [cm]
284 const double rWire = 25.e-4;
285// Outer radius of the tube [cm]
286 const double rTube = 1.451;
287// Half-length of the tube [cm]
288 const double lTube = 10.;
289
290 fGeometrySimple = new Garfield::GeometrySimple();
291// Make a tube (centered at the origin, inner radius: 0, outer radius: rTube).
292 fTube = new Garfield::SolidTube(0., 0., 0, rWire, rTube, lTube);
293// Add the solid to the geometry, together with the medium inside.
294 fGeometrySimple->AddSolid(fTube, fMediumMagboltz);
295 fComponentAnalyticField->SetGeometry(fGeometrySimple);
296
297// Voltages
298 const double vWire = 1000.;
299 const double vTube = 0.;
300// Add the wire in the center.
301 fComponentAnalyticField->AddWire(0., 0., 2 * rWire, vWire, "w");
302// Add the tube.
303 fComponentAnalyticField->AddTube(rTube, vTube, 0, "t");
304
305 fSensor->AddComponent(fComponentAnalyticField);
306
307}
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)
void AddTube(const double radius, const double voltage, const int nEdges, const std::string label)
virtual void SetGeometry(GeometryBase *geo)
void AddSolid(Solid *s, Medium *m)
void AddComponent(ComponentBase *comp)
Definition: Sensor.cc:302

Referenced by InitializePhysics().

◆ DeleteSecondaryParticles()

void GarfieldPhysics::DeleteSecondaryParticles ( )

Definition at line 466 of file GarfieldPhysics.cc.

466 {
467 if (!fSecondaryParticles->empty()) {
468 fSecondaryParticles->erase(fSecondaryParticles->begin(),
469 fSecondaryParticles->end());
470 }
471}

Referenced by DoIt().

◆ Dispose()

void GarfieldPhysics::Dispose ( )
static

Definition at line 47 of file GarfieldPhysics.cc.

47 {
48 delete fGarfieldPhysics;
49 fGarfieldPhysics = 0;
50}

Referenced by main().

◆ DoIt()

void GarfieldPhysics::DoIt ( std::string  particleName,
double  ekin_MeV,
double  time,
double  x_cm,
double  y_cm,
double  z_cm,
double  dx,
double  dy,
double  dz 
)

Definition at line 309 of file GarfieldPhysics.cc.

311 {
312
313 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
314 fEnergyDeposit = 0;
316
317// Wire radius [cm]
318 const double rWire = 25.e-4;
319// Outer radius of the tube [cm]
320 const double rTube = 1.45;
321// Half-length of the tube [cm]
322 const double lTube = 10.;
323
324 double eKin_eV = ekin_MeV * 1e+6;
325
326 double xc = 0., yc = 0., zc = 0., tc = 0.;
327// Number of electrons produced in a collision
328 int nc = 0;
329// Energy loss in a collision
330 double ec = 0.;
331// Dummy variable (not used at present)
332 double extra = 0.;
333 fEnergyDeposit = 0;
334 if (fIonizationModel != "Heed" || particleName == "gamma") {
335 if (particleName == "gamma") {
336 fTrackHeed->TransportPhoton(x_cm, y_cm, z_cm, time, eKin_eV, dx, dy,
337 dz, nc);
338 } else {
339 fTrackHeed->TransportDeltaElectron(x_cm, y_cm, z_cm, time, eKin_eV,
340 dx, dy, dz, nc);
341 fEnergyDeposit = eKin_eV;
342 }
343
344 for (int cl = 0; cl < nc; cl++) {
345 double xe, ye, ze, te;
346 double ee, dxe, dye, dze;
347 fTrackHeed->GetElectron(cl, xe, ye, ze, te, ee, dxe, dye, dze);
348 if (ze < lTube && ze > -lTube && sqrt(xe * xe + ye * ye) < rTube) {
349 nsum++;
350 if (particleName == "gamma") {
351 fEnergyDeposit += fTrackHeed->GetW();
352 }
353 analysisManager->FillH3(1, ze * 10, xe * 10, ye * 10);
354 if (createSecondariesInGeant4) {
355 double newTime = te;
356 if (newTime < time) {
357 newTime += time;
358 }
359 fSecondaryParticles->push_back(
360 new GarfieldParticle("e-",ee, newTime, xe, ye, ze, dxe,
361 dye, dze));
362 }
363
364 fDrift->DriftElectron(xe, ye, ze, te);
365
366 double xe1, ye1, ze1, te1;
367 double xe2, ye2, ze2, te2;
368
369 int status;
370 fDrift->GetElectronEndpoint(0, xe1, ye1, ze1, te1, xe2, ye2,
371 ze2, te2, status);
372
373 if (0 < xe2 && xe2 < rWire) {
374 xe2 += 2 * rWire;
375 }
376 if (0 > xe2 && xe2 > -rWire) {
377 xe2 += -2 * rWire;
378 }
379 if (0 < ye2 && ye2 < rWire) {
380 ye2 += 2 * rWire;
381 }
382 if (0 > ye2 && ye2 > -rWire) {
383 ye2 += -2 * rWire;
384 }
385
386 double e2 = 0.1;
387 fAvalanche->AvalancheElectron(xe2, ye2, ze2, te2, e2, 0, 0, 0);
388
389 int ne = 0, ni = 0;
390 fAvalanche->GetAvalancheSize(ne, ni);
391 fAvalancheSize += ne;
392
393 }
394 }
395 } else {
396 fTrackHeed->SetParticle(particleName);
397 fTrackHeed->SetKineticEnergy(eKin_eV);
398 fTrackHeed->NewTrack(x_cm, y_cm, z_cm, time, dx, dy, dz);
399
400 while (fTrackHeed->GetCluster(xc, yc, zc, tc, nc, ec, extra)) {
401 if (zc < lTube && zc > -lTube && sqrt(xc * xc + yc * yc) < rTube) {
402 nsum += nc;
403 fEnergyDeposit += ec;
404 for (int cl = 0; cl < nc; cl++) {
405 double xe, ye, ze, te;
406 double ee, dxe, dye, dze;
407 fTrackHeed->GetElectron(cl, xe, ye, ze, te, ee, dxe, dye,
408 dze);
409 if (ze < lTube && ze > -lTube
410 && sqrt(xe * xe + ye * ye) < rTube) {
411 analysisManager->FillH3(1, ze * 10, xe * 10, ye * 10);
412 if (createSecondariesInGeant4) {
413 double newTime = te;
414 if (newTime < time) {
415 newTime += time;
416 }
417 fSecondaryParticles->push_back(
418 new GarfieldParticle("e-", ee, newTime, xe, ye,
419 ze, dxe, dye, dze));
420 }
421
422 fDrift->DriftElectron(xe, ye, ze, te);
423
424 double xe1, ye1, ze1, te1;
425 double xe2, ye2, ze2, te2;
426
427 int status;
428 fDrift->GetElectronEndpoint(0, xe1, ye1, ze1, te1, xe2,
429 ye2, ze2, te2, status);
430
431 if (0 < xe2 && xe2 < rWire) {
432 xe2 += 2 * rWire;
433 }
434 if (0 > xe2 && xe2 > -rWire) {
435 xe2 += -2 * rWire;
436 }
437 if (0 < ye2 && ye2 < rWire) {
438 ye2 += 2 * rWire;
439 }
440 if (0 > ye2 && ye2 > -rWire) {
441 ye2 += -2 * rWire;
442 }
443
444 double e2 = 0.1;
445 fAvalanche->AvalancheElectron(xe2, ye2, ze2, te2, e2, 0,
446 0, 0);
447
448 int ne = 0, ni = 0;
449 fAvalanche->GetAvalancheSize(ne, ni);
450 fAvalancheSize += ne;
451
452 }
453 }
454
455 }
456 }
457 }
458 fGain = fAvalancheSize / nsum;
459
460}
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:313
void DeleteSecondaryParticles()
bool DriftElectron(const double x0, const double y0, const double z0, const double t0)
Definition: AvalancheMC.cc:229
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:206
void GetAvalancheSize(int &ne, int &ni) const
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.)
bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra)
Definition: TrackHeed.cc:383
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 &nel)
Definition: TrackHeed.cc:615
bool GetElectron(const int i, double &x, double &y, double &z, double &t, double &e, double &dx, double &dy, double &dz)
Definition: TrackHeed.cc:566
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 &nel)
Definition: TrackHeed.cc:743
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0)
Definition: TrackHeed.cc:158
double GetW() const
Definition: TrackHeed.cc:1370
void SetKineticEnergy(const double ekin)
Definition: Track.cc:171
virtual void SetParticle(std::string part)
Definition: Track.cc:29

Referenced by GarfieldG4FastSimulationModel::DoIt().

◆ EnableCreateSecondariesInGeant4()

void GarfieldPhysics::EnableCreateSecondariesInGeant4 ( bool  flag)
inline

Definition at line 96 of file GarfieldPhysics.hh.

96{createSecondariesInGeant4 = flag;};

◆ FindParticleName()

bool GarfieldPhysics::FindParticleName ( const std::string  name,
std::string  program = "garfield" 
)

Definition at line 168 of file GarfieldPhysics.cc.

168 {
169 MapParticlesEnergy::iterator it;
170 if (program == "garfield") {
171 it = fMapParticlesEnergyGarfield->find(name);
172 if (it != fMapParticlesEnergyGarfield->end()) {
173 return true;
174 }
175 return false;
176 } else {
177 it = fMapParticlesEnergyGeant4->find(name);
178 if (it != fMapParticlesEnergyGeant4->end()) {
179 return true;
180 }
181 return false;
182 }
183}

Referenced by GarfieldPhysicsList::AddParameterisation(), and GarfieldG4FastSimulationModel::IsApplicable().

◆ FindParticleNameEnergy()

bool GarfieldPhysics::FindParticleNameEnergy ( std::string  name,
double  ekin_MeV,
std::string  program = "garfield" 
)

Definition at line 185 of file GarfieldPhysics.cc.

186 {
187 MapParticlesEnergy::iterator it;
188 if (program == "garfield") {
189 it = fMapParticlesEnergyGarfield->find(name);
190 if (it != fMapParticlesEnergyGarfield->end()) {
191 EnergyRange_MeV range = it->second;
192 if (range.first <= ekin_MeV && range.second >= ekin_MeV) {
193 return true;
194 }
195 }
196 return false;
197 } else {
198 it = fMapParticlesEnergyGeant4->find(name);
199 if (it != fMapParticlesEnergyGeant4->end()) {
200 EnergyRange_MeV range = it->second;
201 if (range.first <= ekin_MeV && range.second >= ekin_MeV) {
202 return true;
203 }
204 }
205 return false;
206 }
207}
std::pair< double, double > EnergyRange_MeV

Referenced by GarfieldG4FastSimulationModel::ModelTrigger().

◆ GetAvalancheSize()

double GarfieldPhysics::GetAvalancheSize ( )
inline

Definition at line 99 of file GarfieldPhysics.hh.

99{return fAvalancheSize;};

Referenced by GarfieldEventAction::EndOfEventAction().

◆ GetCreateSecondariesInGeant4()

bool GarfieldPhysics::GetCreateSecondariesInGeant4 ( )
inline

Definition at line 97 of file GarfieldPhysics.hh.

97{return createSecondariesInGeant4;};

Referenced by GarfieldG4FastSimulationModel::DoIt().

◆ GetEnergyDeposit_MeV()

double GarfieldPhysics::GetEnergyDeposit_MeV ( )
inline

Definition at line 98 of file GarfieldPhysics.hh.

98{return fEnergyDeposit/1000000;};

Referenced by GarfieldG4FastSimulationModel::DoIt().

◆ GetGain()

double GarfieldPhysics::GetGain ( )
inline

Definition at line 100 of file GarfieldPhysics.hh.

100{return fGain;};

Referenced by GarfieldEventAction::EndOfEventAction().

◆ GetInstance()

GarfieldPhysics * GarfieldPhysics::GetInstance ( )
static

◆ GetIonizationModel()

std::string GarfieldPhysics::GetIonizationModel ( )

Definition at line 87 of file GarfieldPhysics.cc.

87 {
88 return fIonizationModel;
89}

Referenced by GarfieldPhysicsList::AddParameterisation().

◆ GetMaxEnergyMeVParticle()

double GarfieldPhysics::GetMaxEnergyMeVParticle ( std::string  name,
std::string  program = "garfield" 
)

Definition at line 229 of file GarfieldPhysics.cc.

230 {
231 MapParticlesEnergy::iterator it;
232 if (program == "garfield") {
233 it = fMapParticlesEnergyGarfield->find(name);
234 if (it != fMapParticlesEnergyGarfield->end()) {
235 EnergyRange_MeV range = it->second;
236 return range.second;
237
238 }
239 } else {
240 it = fMapParticlesEnergyGeant4->find(name);
241 if (it != fMapParticlesEnergyGeant4->end()) {
242 EnergyRange_MeV range = it->second;
243 return range.second;
244 }
245 }
246 return -1;
247}

Referenced by GarfieldPhysicsList::AddParameterisation().

◆ GetMinEnergyMeVParticle()

double GarfieldPhysics::GetMinEnergyMeVParticle ( std::string  name,
std::string  program = "garfield" 
)

Definition at line 209 of file GarfieldPhysics.cc.

210 {
211 MapParticlesEnergy::iterator it;
212 if (program == "garfield") {
213 it = fMapParticlesEnergyGarfield->find(name);
214 if (it != fMapParticlesEnergyGarfield->end()) {
215 EnergyRange_MeV range = it->second;
216 return range.first;
217 }
218 return false;
219 } else {
220 it = fMapParticlesEnergyGeant4->find(name);
221 if (it != fMapParticlesEnergyGeant4->end()) {
222 EnergyRange_MeV range = it->second;
223 return range.first;
224 }
225 }
226 return -1;
227}

Referenced by GarfieldPhysicsList::AddParameterisation().

◆ GetSecondaryParticles()

std::vector< GarfieldParticle * > * GarfieldPhysics::GetSecondaryParticles ( )

Definition at line 462 of file GarfieldPhysics.cc.

462 {
463 return fSecondaryParticles;
464}

Referenced by GarfieldG4FastSimulationModel::DoIt().

◆ InitializePhysics()

void GarfieldPhysics::InitializePhysics ( )

Definition at line 249 of file GarfieldPhysics.cc.

249 {
250
251 fMediumMagboltz = new Garfield::MediumMagboltz();
252
253 fMediumMagboltz->SetComposition("ar", 70., "co2", 30.);
254 fMediumMagboltz->SetTemperature(293.15);
255 fMediumMagboltz->SetPressure(760.);
256 fMediumMagboltz->EnableDebugging();
257 fMediumMagboltz->Initialise(true);
258 fMediumMagboltz->DisableDebugging();
259// Set the Penning transfer efficiency.
260 const double rPenning = 0.57;
261 const double lambdaPenning = 0.;
262 fMediumMagboltz->EnablePenningTransfer(rPenning, lambdaPenning, "ar");
263 fMediumMagboltz->LoadGasFile("ar_70_co2_30_1000mbar.gas");
264
265 fSensor = new Garfield::Sensor();
266 fDrift = new Garfield::AvalancheMC();
267 fAvalanche = new Garfield::AvalancheMicroscopic();
268 fComponentAnalyticField = new Garfield::ComponentAnalyticField();
269
271
272 fDrift->SetSensor(fSensor);
273 fAvalanche->SetSensor(fSensor);
274
275 fTrackHeed = new Garfield::TrackHeed();
276 fTrackHeed->EnableDebugging();
277 fTrackHeed->SetSensor(fSensor);
278
279 fTrackHeed->EnableDeltaElectronTransport();
280}
void SetSensor(Sensor *s)
Definition: AvalancheMC.cc:52
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.)
Definition: MediumGas.cc:64
bool LoadGasFile(const std::string &filename)
Definition: MediumGas.cc:283
void EnablePenningTransfer(const double r, const double lambda)
bool Initialise(const bool verbose=false)
void SetTemperature(const double &t)
Definition: Medium.cc:117
void SetPressure(const double &p)
Definition: Medium.cc:128
void DisableDebugging()
Definition: Medium.hh:281
void EnableDebugging()
Definition: Medium.hh:280
void EnableDeltaElectronTransport()
Definition: TrackHeed.hh:67
void SetSensor(Sensor *s)
Definition: Track.cc:185
void EnableDebugging()
Definition: Track.hh:57

Referenced by GarfieldG4FastSimulationModel::GarfieldG4FastSimulationModel().

◆ SetIonizationModel()

void GarfieldPhysics::SetIonizationModel ( std::string  model,
bool  useDefaults = true 
)

Definition at line 91 of file GarfieldPhysics.cc.

91 {
92 if (model != "PAIPhot" && model != "PAI" && model != "Heed") {
93
94 std::cout << "Unknown ionization model " << model << std::endl;
95 std::cout << "Using PAIPhot as default model!" << std::endl;
96 model = "PAIPhot";
97 }
98 fIonizationModel = model;
99
100 if (fIonizationModel == "PAIPhot" || fIonizationModel == "PAI") {
101 if (useDefaults == true) {
102 //Particle types and energies for which the G4FastSimulationModel with Garfield++ is valid
103 this->AddParticleName("e-", 1e-6, 1e-3, "garfield");
104 this->AddParticleName("gamma", 1e-6, 1e+8, "garfield");
105
106 //Particle types and energies for which the PAI or PAIPhot model is valid
107 this->AddParticleName("e-", 0, 1e+8, "geant4");
108 this->AddParticleName("e+", 0, 1e+8, "geant4");
109 this->AddParticleName("mu-", 0, 1e+8, "geant4");
110 this->AddParticleName("mu+", 0, 1e+8, "geant4");
111 this->AddParticleName("proton", 0, 1e+8, "geant4");
112 this->AddParticleName("pi+", 0, 1e+8, "geant4");
113 this->AddParticleName("pi-", 0, 1e+8, "geant4");
114 this->AddParticleName("alpha", 0, 1e+8, "geant4");
115 this->AddParticleName("He3", 0, 1e+8, "geant4");
116 this->AddParticleName("GenericIon", 0, 1e+8, "geant4");
117 }
118
119 } else if (fIonizationModel == "Heed") {
120 if (useDefaults == true) {
121 //Particle types and energies for which the G4FastSimulationModel with Garfield++ is valid
122 this->AddParticleName("gamma", 1e-6, 1e+8, "garfield");
123 this->AddParticleName("e-", 6e-2, 1e+7, "garfield");
124 this->AddParticleName("e+", 6e-2, 1e+7, "garfield");
125 this->AddParticleName("mu-", 1e+1, 1e+8, "garfield");
126 this->AddParticleName("mu+", 1e+1, 1e+8, "garfield");
127 this->AddParticleName("pi-", 2e+1, 1e+8, "garfield");
128 this->AddParticleName("pi+", 2e+1, 1e+8, "garfield");
129 this->AddParticleName("kaon-", 1e+1, 1e+8, "garfield");
130 this->AddParticleName("kaon+", 1e+1, 1e+8, "garfield");
131 this->AddParticleName("proton", 9.e+1, 1e+8, "garfield");
132 this->AddParticleName("anti_proton", 9.e+1, 1e+8, "garfield");
133 this->AddParticleName("deuteron", 2.e+2, 1e+8, "garfield");
134 this->AddParticleName("alpha", 4.e+2, 1e+8, "garfield");
135 }
136
137 }
138}
void AddParticleName(const std::string particleName, double ekin_min_MeV, double ekin_max_MeV, std::string program)

Referenced by GarfieldMessenger::SetNewValue().


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