32#include "TGeoManager.h"
41 if (!fGarfieldPhysics) {
44 return fGarfieldPhysics;
48 delete fGarfieldPhysics;
52GarfieldPhysics::GarfieldPhysics() {
55 fSecondaryParticles =
new std::vector<GarfieldParticle*>();
60 fComponentAnalyticField = 0;
66 createSecondariesInGeant4 =
false;
67 fIonizationModel =
"PAIPhot";
70GarfieldPhysics::~GarfieldPhysics() {
71 delete fMapParticlesEnergyGeant4;
72 delete fMapParticlesEnergyGarfield;
74 delete fSecondaryParticles;
75 delete fMediumMagboltz;
79 delete fComponentAnalyticField;
82 delete fGeometrySimple;
84 std::cout <<
"Deconstructor GarfieldPhysics" << std::endl;
88 return fIonizationModel;
92 if (model !=
"PAIPhot" && model !=
"PAI" && model !=
"Heed") {
94 std::cout <<
"Unknown ionization model " << model << std::endl;
95 std::cout <<
"Using PAIPhot as default model!" << std::endl;
98 fIonizationModel = model;
100 if (fIonizationModel ==
"PAIPhot" || fIonizationModel ==
"PAI") {
101 if (useDefaults ==
true) {
119 }
else if (fIonizationModel ==
"Heed") {
120 if (useDefaults ==
true) {
141 double ekin_min_MeV,
double ekin_max_MeV, std::string program) {
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"
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;
154 fMapParticlesEnergyGarfield->insert(
155 std::make_pair(particleName,
156 std::make_pair(ekin_min_MeV, ekin_max_MeV)));
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)));
169 MapParticlesEnergy::iterator it;
170 if (program ==
"garfield") {
171 it = fMapParticlesEnergyGarfield->find(name);
172 if (it != fMapParticlesEnergyGarfield->end()) {
177 it = fMapParticlesEnergyGeant4->find(name);
178 if (it != fMapParticlesEnergyGeant4->end()) {
186 std::string program) {
187 MapParticlesEnergy::iterator it;
188 if (program ==
"garfield") {
189 it = fMapParticlesEnergyGarfield->find(name);
190 if (it != fMapParticlesEnergyGarfield->end()) {
192 if (range.first <= ekin_MeV && range.second >= ekin_MeV) {
198 it = fMapParticlesEnergyGeant4->find(name);
199 if (it != fMapParticlesEnergyGeant4->end()) {
201 if (range.first <= ekin_MeV && range.second >= ekin_MeV) {
210 std::string program) {
211 MapParticlesEnergy::iterator it;
212 if (program ==
"garfield") {
213 it = fMapParticlesEnergyGarfield->find(name);
214 if (it != fMapParticlesEnergyGarfield->end()) {
220 it = fMapParticlesEnergyGeant4->find(name);
221 if (it != fMapParticlesEnergyGeant4->end()) {
230 std::string program) {
231 MapParticlesEnergy::iterator it;
232 if (program ==
"garfield") {
233 it = fMapParticlesEnergyGarfield->find(name);
234 if (it != fMapParticlesEnergyGarfield->end()) {
240 it = fMapParticlesEnergyGeant4->find(name);
241 if (it != fMapParticlesEnergyGeant4->end()) {
260 const double rPenning = 0.57;
261 const double lambdaPenning = 0.;
263 fMediumMagboltz->
LoadGasFile(
"ar_70_co2_30_1000mbar.gas");
284 const double rWire = 25.e-4;
286 const double rTube = 1.451;
288 const double lTube = 10.;
294 fGeometrySimple->
AddSolid(fTube, fMediumMagboltz);
295 fComponentAnalyticField->
SetGeometry(fGeometrySimple);
298 const double vWire = 1000.;
299 const double vTube = 0.;
301 fComponentAnalyticField->
AddWire(0., 0., 2 * rWire, vWire,
"w");
303 fComponentAnalyticField->
AddTube(rTube, vTube, 0,
"t");
310 double time,
double x_cm,
double y_cm,
double z_cm,
double dx,
311 double dy,
double dz) {
313 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
318 const double rWire = 25.e-4;
320 const double rTube = 1.45;
322 const double lTube = 10.;
324 double eKin_eV = ekin_MeV * 1e+6;
326 double xc = 0., yc = 0., zc = 0., tc = 0.;
334 if (fIonizationModel !=
"Heed" || particleName ==
"gamma") {
335 if (particleName ==
"gamma") {
341 fEnergyDeposit = eKin_eV;
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) {
350 if (particleName ==
"gamma") {
351 fEnergyDeposit += fTrackHeed->
GetW();
353 analysisManager->FillH3(1, ze * 10, xe * 10, ye * 10);
354 if (createSecondariesInGeant4) {
356 if (newTime < time) {
359 fSecondaryParticles->push_back(
366 double xe1, ye1, ze1, te1;
367 double xe2, ye2, ze2, te2;
373 if (0 < xe2 && xe2 < rWire) {
376 if (0 > xe2 && xe2 > -rWire) {
379 if (0 < ye2 && ye2 < rWire) {
382 if (0 > ye2 && ye2 > -rWire) {
391 fAvalancheSize += ne;
398 fTrackHeed->
NewTrack(x_cm, y_cm, z_cm, time, dx, dy, dz);
400 while (fTrackHeed->
GetCluster(xc, yc, zc, tc, nc, ec, extra)) {
401 if (zc < lTube && zc > -lTube &&
sqrt(xc * xc + yc * yc) < rTube) {
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,
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) {
414 if (newTime < time) {
417 fSecondaryParticles->push_back(
424 double xe1, ye1, ze1, te1;
425 double xe2, ye2, ze2, te2;
429 ye2, ze2, te2, status);
431 if (0 < xe2 && xe2 < rWire) {
434 if (0 > xe2 && xe2 > -rWire) {
437 if (0 < ye2 && ye2 < rWire) {
440 if (0 > ye2 && ye2 > -rWire) {
450 fAvalancheSize += ne;
458 fGain = fAvalancheSize / nsum;
463 return fSecondaryParticles;
467 if (!fSecondaryParticles->empty()) {
468 fSecondaryParticles->erase(fSecondaryParticles->begin(),
469 fSecondaryParticles->end());
DoubleAc sqrt(const DoubleAc &f)
Selection of the analysis technology.
Definition of the GarfieldPhysics class.
std::map< const std::string, EnergyRange_MeV > MapParticlesEnergy
std::pair< double, double > EnergyRange_MeV
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)
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
void SetSensor(Sensor *s)
void GetAvalancheSize(int &ne, int &ni) const
void SetSensor(Sensor *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.)
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)
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.)
bool LoadGasFile(const std::string &filename)
void EnablePenningTransfer(const double r, const double lambda)
bool Initialise(const bool verbose=false)
void SetTemperature(const double &t)
void SetPressure(const double &p)
void AddComponent(ComponentBase *comp)
bool GetCluster(double &xcls, double &ycls, double &zcls, double &tcls, int &n, double &e, double &extra)
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)
bool GetElectron(const int i, double &x, double &y, double &z, double &t, double &e, double &dx, double &dy, double &dz)
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)
void EnableDeltaElectronTransport()
bool NewTrack(const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0)
void SetSensor(Sensor *s)
void SetKineticEnergy(const double ekin)
virtual void SetParticle(std::string part)