31#ifndef GarfieldPhysics_h
32#define GarfieldPhysics_h 1
54 GarfieldParticle(std::string particleName,
double ekin_eV,
double time,
double x_cm,
double y_cm,
double z_cm,
double dx,
double dy,
double dz):fParticleName(particleName), fEkin_MeV(ekin_eV/1000000), fTime(time), fx_mm(10*x_cm),fy_mm(10*y_cm), fz_mm(10*z_cm), fdx(dx), fdy(dy), fdz(dz){}
69 std::string fParticleName;
70 double fEkin_MeV, fTime, fx_mm,fy_mm,fz_mm,fdx,fdy,fdz;
84 void DoIt(std::string particleName,
double ekin_MeV,
double time,
85 double x_cm,
double y_cm,
double z_cm,
double dx,
double dy,
double dz);
87 void AddParticleName(
const std::string particleName,
double ekin_min_MeV,
double ekin_max_MeV, std::string program);
88 bool FindParticleName(
const std::string name, std::string program =
"garfield");
101 inline void Clear() {fEnergyDeposit=0;fAvalancheSize=0;fGain=0;nsum=0;}
107 std::string fIonizationModel;
112 TGeoManager* fGeoManager;
123 std::vector<GarfieldParticle*>* fSecondaryParticles;
125 bool createSecondariesInGeant4;
126 double fEnergyDeposit;
127 double fAvalancheSize;
std::map< const std::string, EnergyRange_MeV > MapParticlesEnergy
std::pair< double, double > EnergyRange_MeV
std::string getParticleName()
GarfieldParticle(std::string particleName, double ekin_eV, double time, double x_cm, double y_cm, double z_cm, double dx, double dy, double dz)
double GetMaxEnergyMeVParticle(std::string name, std::string program="garfield")
bool GetCreateSecondariesInGeant4()
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()
double GetEnergyDeposit_MeV()
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")
void EnableCreateSecondariesInGeant4(bool flag)
double GetAvalancheSize()
bool FindParticleName(const std::string name, std::string program="garfield")
Calculate electron drift lines and avalanches using microscopic tracking.
Use a geometry defined using the ROOT TGeo package.
"Native" geometry, using simple shapes.
Generate tracks using Heed++.