Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::TrackHeed Class Reference

Generate tracks using Heed++. More...

#include <TrackHeed.hh>

Inheritance diagram for Garfield::TrackHeed:
Garfield::Track

Classes

struct  Cluster
struct  Electron
struct  Photon

Public Member Functions

 TrackHeed ()
 Default constructor.
 TrackHeed (Sensor *sensor)
 Constructor.
virtual ~TrackHeed ()
 Destructor.
bool NewTrack (const double x0, const double y0, const double z0, const double t0, const double dx0, const double dy0, const double dz0) override
 Calculate a new track starting from (x0, y0, z0) at time t0 in direction (dx0, dy0, dz0).
const std::vector< Cluster > & GetClusters () const
bool GetCluster (double &xc, double &yc, double &zc, double &tc, int &nc, double &ec, double &extra)
bool GetCluster (double &xc, double &yc, double &zc, double &tc, int &ne, int &ni, double &ec, double &extra)
bool GetCluster (double &xc, double &yc, double &zc, double &tc, int &ne, int &ni, int &np, double &ec, double &extra)
 Get the next "cluster" (ionising collision of the charged particle).
bool GetElectron (const unsigned int i, double &x, double &y, double &z, double &t, double &e, double &dx, double &dy, double &dz)
 Retrieve the properties of a conduction or delta electron in the current cluster.
bool GetIon (const unsigned int i, double &x, double &y, double &z, double &t) const
 Retrieve the properties of an ion in the current cluster.
bool GetPhoton (const unsigned int i, double &x, double &y, double &z, double &t, double &e, double &dx, double &dy, double &dz) const
 Retrieve the properties of an unabsorbed photon.
double GetClusterDensity () override
 Get the cluster density (number of ionizing collisions per cm or inverse mean free path for ionization).
double GetStoppingPower () override
 Get the stopping power (mean energy loss [eV] per cm).
double GetW () const
 Return the W value of the medium (of the last simulated track).
double GetFanoFactor () const
 Return the Fano factor of the medium (of the last simulated track).
double GetPhotoAbsorptionCrossSection (const double e) const
 Return the photoabsorption cross-section at a given energy.
bool Initialise (Medium *medium, const bool verbose=false)
 Compute the differential cross-section for a given medium.
Cluster 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)
 Simulate a delta electron.
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)
 Simulate a delta electron.
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)
 Simulate a delta electron.
Cluster 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)
 Simulate a photon.
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, int &np)
 Simulate a photon.
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)
 Simulate a photon.
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)
 Simulate a photon.
void EnableElectricField ()
 Take the electric field into account in the stepping algorithm.
void DisableElectricField ()
 Do not take the electric field into account in the stepping algorithm.
void EnableMagneticField ()
 Take the magnetic field into account in the stepping algorithm.
void DisableMagneticField ()
 Do not take the magnetic field into account in the stepping algorithm.
void SetSteppingLimits (const double maxStep, const double radStraight, const double stepAngleStraight, const double stepAngleCurved)
 Set parameters for calculating the particle trajectory.
void GetSteppingLimits (double &maxStep, double &radStraight, double &stepAngleStraight, double &stepAngleCurved)
void CrossInactiveMedia (const bool on=true)
void EnableCoulombScattering (const bool on=true)
void EnableDeltaElectronTransport ()
 Switch simulation of delta electrons on.
void DisableDeltaElectronTransport ()
 Switch simulation of delta electrons off.
void EnablePhotonReabsorption (const bool on=true)
 Simulate (or not) the photons produced in the atomic relaxation cascade.
void EnablePhotoAbsorptionCrossSectionOutput (const bool on)
 Write the photoabsorption cross-sections used to a text file.
void SetEnergyMesh (const double e0, const double e1, const int nsteps)
 Specify the energy mesh to be used.
void SetParticleUser (const double m, const double z)
 Define particle mass and charge (for exotic particles).
void EnableOneStepFly (const bool on)
Public Member Functions inherited from Garfield::Track
 Track ()=delete
 Default constructor.
 Track (const std::string &name)
 Constructor.
virtual ~Track ()
 Destructor.
virtual void SetParticle (const std::string &part)
 Set the type of charged particle.
void SetEnergy (const double e)
 Set the particle energy.
void SetBetaGamma (const double bg)
 Set the relative momentum of the particle.
void SetBeta (const double beta)
 Set the speed ( $\beta = v/c$) of the particle.
void SetGamma (const double gamma)
 Set the Lorentz factor of the particle.
void SetMomentum (const double p)
 Set the particle momentum.
void SetKineticEnergy (const double ekin)
 Set the kinetic energy of the particle.
double GetEnergy () const
 Return the particle energy.
double GetBetaGamma () const
 Return the $\beta\gamma$ of the projectile.
double GetBeta () const
 Return the speed ( $\beta = v/c$) of the projectile.
double GetGamma () const
 Return the Lorentz factor of the projectile.
double GetMomentum () const
 Return the particle momentum.
double GetKineticEnergy () const
 Return the kinetic energy of the projectile.
double GetCharge () const
 Get the charge of the projectile.
double GetMass () const
 Get the mass [eV / c2] of the projectile.
void SetSensor (Sensor *s)
 Set the sensor through which to transport the particle.
void EnablePlotting (ViewDrift *viewer)
 Switch on plotting.
void DisablePlotting ()
 Switch off plotting.
void EnableDebugging ()
 Switch on debugging messages.
void DisableDebugging ()
 Switch off debugging messages.

Private Member Functions

 TrackHeed (const TrackHeed &heed)
TrackHeedoperator= (const TrackHeed &heed)
bool SetupGas (Medium *medium)
bool SetupMaterial (Medium *medium)
bool SetupDelta (const std::string &databasePath)
bool AddCluster (Heed::HeedPhoton *virtualPhoton, std::vector< Cluster > &clusters)
void AddElectrons (const std::vector< Heed::HeedCondElectron > &conductionElectrons, std::vector< Electron > &electrons)
bool IsInside (const double x, const double y, const double z)
bool UpdateBoundingBox (bool &update)

Private Attributes

bool m_oneStepFly = false
bool m_ready = false
bool m_hasActiveTrack = false
double m_mediumDensity = -1.
std::string m_mediumName = ""
bool m_usePacsOutput = false
bool m_doPhotonReabsorption = true
bool m_crossInactiveMedia = false
bool m_coulombScattering = false
bool m_useBfieldAuto = true
bool m_doDeltaTransport = true
std::vector< Clusterm_clusters
size_t m_cluster = 0
std::unique_ptr< Heed::particle_def > m_particle_def
std::unique_ptr< Heed::HeedMatterDef > m_matter
std::unique_ptr< Heed::GasDef > m_gas
std::unique_ptr< Heed::MatterDef > m_material
double m_emin = 2.e-6
double m_emax = 2.e-1
unsigned int m_nEnergyIntervals = 200
std::unique_ptr< Heed::EnergyMesh > m_energyMesh
std::unique_ptr< Heed::EnTransfCS > m_transferCs
std::unique_ptr< Heed::ElElasticScat > m_elScat
std::unique_ptr< Heed::ElElasticScatLowSigma > m_lowSigma
std::unique_ptr< Heed::PairProd > m_pairProd
std::unique_ptr< Heed::HeedDeltaElectronCS > m_deltaCs
std::unique_ptr< HeedChamber > m_chamber
std::unique_ptr< HeedFieldMap > m_fieldMap
double m_lX = 0.
double m_lY = 0.
double m_lZ = 0.
double m_cX = 0.
double m_cY = 0.
double m_cZ = 0.
double m_maxStep = 100.
 Max. step length.
double m_radStraight = 1000.
 Bending radius beyond which to use straight-line approximation.
double m_stepAngleStraight = 0.1
 Angular step for curved trajectories approximated by straight-line steps.
double m_stepAngleCurved = 0.2
 Angular step for curved lines.

Additional Inherited Members

Protected Member Functions inherited from Garfield::Track
void PlotNewTrack (const double x0, const double y0, const double z0)
void PlotCluster (const double x0, const double y0, const double z0)
Static Protected Member Functions inherited from Garfield::Track
static std::array< double, 3 > StepBfield (const double dt, const double qoverm, const double vmag, double bx, double by, double bz, std::array< double, 3 > &dir)
Protected Attributes inherited from Garfield::Track
std::string m_className = "Track"
double m_q = -1.
int m_spin = 1
double m_mass
double m_energy = 0.
double m_beta2 = 1.
bool m_isElectron = false
std::string m_particleName = "mu-"
Sensor * m_sensor = nullptr
bool m_isChanged = true
ViewDriftm_viewer = nullptr
bool m_debug = false
size_t m_plotId = 0

Detailed Description

Generate tracks using Heed++.

Definition at line 34 of file TrackHeed.hh.

Constructor & Destructor Documentation

◆ TrackHeed() [1/3]

Garfield::TrackHeed::TrackHeed ( )
inline

Default constructor.

Definition at line 68 of file TrackHeed.hh.

68: TrackHeed(nullptr) {}
TrackHeed()
Default constructor.
Definition TrackHeed.hh:68

◆ TrackHeed() [2/3]

Garfield::TrackHeed::TrackHeed ( Sensor * sensor)

Constructor.

◆ ~TrackHeed()

virtual Garfield::TrackHeed::~TrackHeed ( )
virtual

Destructor.

◆ TrackHeed() [3/3]

Garfield::TrackHeed::TrackHeed ( const TrackHeed & heed)
private

Member Function Documentation

◆ AddCluster()

bool Garfield::TrackHeed::AddCluster ( Heed::HeedPhoton * virtualPhoton,
std::vector< Cluster > & clusters )
private

◆ AddElectrons()

void Garfield::TrackHeed::AddElectrons ( const std::vector< Heed::HeedCondElectron > & conductionElectrons,
std::vector< Electron > & electrons )
private

◆ CrossInactiveMedia()

void Garfield::TrackHeed::CrossInactiveMedia ( const bool on = true)
inline

Definition at line 250 of file TrackHeed.hh.

250{ m_crossInactiveMedia = on; }

◆ DisableDeltaElectronTransport()

void Garfield::TrackHeed::DisableDeltaElectronTransport ( )
inline

Switch simulation of delta electrons off.

Definition at line 257 of file TrackHeed.hh.

257{ m_doDeltaTransport = false; }

◆ DisableElectricField()

void Garfield::TrackHeed::DisableElectricField ( )

Do not take the electric field into account in the stepping algorithm.

◆ DisableMagneticField()

void Garfield::TrackHeed::DisableMagneticField ( )

Do not take the magnetic field into account in the stepping algorithm.

◆ EnableCoulombScattering()

void Garfield::TrackHeed::EnableCoulombScattering ( const bool on = true)
inline

Definition at line 251 of file TrackHeed.hh.

251 {
253 }

◆ EnableDeltaElectronTransport()

void Garfield::TrackHeed::EnableDeltaElectronTransport ( )
inline

Switch simulation of delta electrons on.

Definition at line 255 of file TrackHeed.hh.

255{ m_doDeltaTransport = true; }

◆ EnableElectricField()

void Garfield::TrackHeed::EnableElectricField ( )

Take the electric field into account in the stepping algorithm.

◆ EnableMagneticField()

void Garfield::TrackHeed::EnableMagneticField ( )

Take the magnetic field into account in the stepping algorithm.

◆ EnableOneStepFly()

void Garfield::TrackHeed::EnableOneStepFly ( const bool on)
inline

Definition at line 278 of file TrackHeed.hh.

278{ m_oneStepFly = on; }

◆ EnablePhotoAbsorptionCrossSectionOutput()

void Garfield::TrackHeed::EnablePhotoAbsorptionCrossSectionOutput ( const bool on)
inline

Write the photoabsorption cross-sections used to a text file.

Definition at line 265 of file TrackHeed.hh.

265 {
266 m_usePacsOutput = on;
267 }

◆ EnablePhotonReabsorption()

void Garfield::TrackHeed::EnablePhotonReabsorption ( const bool on = true)
inline

Simulate (or not) the photons produced in the atomic relaxation cascade.

Definition at line 260 of file TrackHeed.hh.

260 {
262 }

◆ GetCluster() [1/3]

bool Garfield::TrackHeed::GetCluster ( double & xc,
double & yc,
double & zc,
double & tc,
int & nc,
double & ec,
double & extra )

◆ GetCluster() [2/3]

bool Garfield::TrackHeed::GetCluster ( double & xc,
double & yc,
double & zc,
double & tc,
int & ne,
int & ni,
double & ec,
double & extra )

◆ GetCluster() [3/3]

bool Garfield::TrackHeed::GetCluster ( double & xc,
double & yc,
double & zc,
double & tc,
int & ne,
int & ni,
int & np,
double & ec,
double & extra )

Get the next "cluster" (ionising collision of the charged particle).

Parameters
xc,yc,zccoordinates of the collision
tctime of the collision
nenumber of electrons
ninumber of ions
npnumber of fluorescence photons
ecdeposited energy
extraadditional information (not always implemented)

◆ GetClusterDensity()

double Garfield::TrackHeed::GetClusterDensity ( )
overridevirtual

Get the cluster density (number of ionizing collisions per cm or inverse mean free path for ionization).

Reimplemented from Garfield::Track.

◆ GetClusters()

const std::vector< Cluster > & Garfield::TrackHeed::GetClusters ( ) const
inline

Definition at line 77 of file TrackHeed.hh.

77{ return m_clusters; }
std::vector< Cluster > m_clusters
Definition TrackHeed.hh:301

◆ GetElectron()

bool Garfield::TrackHeed::GetElectron ( const unsigned int i,
double & x,
double & y,
double & z,
double & t,
double & e,
double & dx,
double & dy,
double & dz )

Retrieve the properties of a conduction or delta electron in the current cluster.

Parameters
iindex of the electron
x,y,zcoordinates of the electron
ttime
ekinetic energy (only meaningful for delta-electrons)
dx,dy,dzdirection vector (only meaningful for delta-electrons)

◆ GetFanoFactor()

double Garfield::TrackHeed::GetFanoFactor ( ) const

Return the Fano factor of the medium (of the last simulated track).

◆ GetIon()

bool Garfield::TrackHeed::GetIon ( const unsigned int i,
double & x,
double & y,
double & z,
double & t ) const

Retrieve the properties of an ion in the current cluster.

Parameters
iindex of the ion
x,y,zcoordinates of the ion
ttime

◆ GetPhotoAbsorptionCrossSection()

double Garfield::TrackHeed::GetPhotoAbsorptionCrossSection ( const double e) const

Return the photoabsorption cross-section at a given energy.

◆ GetPhoton()

bool Garfield::TrackHeed::GetPhoton ( const unsigned int i,
double & x,
double & y,
double & z,
double & t,
double & e,
double & dx,
double & dy,
double & dz ) const

Retrieve the properties of an unabsorbed photon.

Parameters
iindex of the photon
x,y,zcoordinates of the photon
ttime
ephoton energy
dx,dy,dzdirection vector

◆ GetSteppingLimits()

void Garfield::TrackHeed::GetSteppingLimits ( double & maxStep,
double & radStraight,
double & stepAngleStraight,
double & stepAngleCurved )
inline

Definition at line 242 of file TrackHeed.hh.

243 {
244 maxStep = m_maxStep;
245 radStraight = m_radStraight;
246 stepAngleStraight = m_stepAngleStraight;
247 stepAngleCurved = m_stepAngleCurved;
248 }
double m_stepAngleCurved
Angular step for curved lines.
Definition TrackHeed.hh:340
double m_stepAngleStraight
Angular step for curved trajectories approximated by straight-line steps.
Definition TrackHeed.hh:338
double m_maxStep
Max. step length.
Definition TrackHeed.hh:334
double m_radStraight
Bending radius beyond which to use straight-line approximation.
Definition TrackHeed.hh:336

◆ GetStoppingPower()

double Garfield::TrackHeed::GetStoppingPower ( )
overridevirtual

Get the stopping power (mean energy loss [eV] per cm).

Reimplemented from Garfield::Track.

◆ GetW()

double Garfield::TrackHeed::GetW ( ) const

Return the W value of the medium (of the last simulated track).

◆ Initialise()

bool Garfield::TrackHeed::Initialise ( Medium * medium,
const bool verbose = false )

Compute the differential cross-section for a given medium.

◆ IsInside()

bool Garfield::TrackHeed::IsInside ( const double x,
const double y,
const double z )
private

◆ NewTrack()

bool Garfield::TrackHeed::NewTrack ( const double x0,
const double y0,
const double z0,
const double t0,
const double dx0,
const double dy0,
const double dz0 )
overridevirtual

Calculate a new track starting from (x0, y0, z0) at time t0 in direction (dx0, dy0, dz0).

Implements Garfield::Track.

◆ operator=()

TrackHeed & Garfield::TrackHeed::operator= ( const TrackHeed & heed)
private

◆ SetEnergyMesh()

void Garfield::TrackHeed::SetEnergyMesh ( const double e0,
const double e1,
const int nsteps )

Specify the energy mesh to be used.

Parameters
e0,e1lower/higher limit of the energy range [eV]
nstepsnumber of intervals

◆ SetParticleUser()

void Garfield::TrackHeed::SetParticleUser ( const double m,
const double z )

Define particle mass and charge (for exotic particles).

For standard particles Track::SetParticle should be used.

◆ SetSteppingLimits()

void Garfield::TrackHeed::SetSteppingLimits ( const double maxStep,
const double radStraight,
const double stepAngleStraight,
const double stepAngleCurved )
inline

Set parameters for calculating the particle trajectory.

Parameters
maxStepmaximum step length
radStraightradius beyond which to approximate circles by polylines.
stepAngleStraightmax. angular step (in radian) when using polyline steps.
stepAngleCurvedmax. angular step (in radian) when using circular steps.

Definition at line 234 of file TrackHeed.hh.

236 {
237 m_maxStep = maxStep;
238 m_radStraight = radStraight;
239 m_stepAngleStraight = stepAngleStraight;
240 m_stepAngleCurved = stepAngleCurved;
241 }

◆ SetupDelta()

bool Garfield::TrackHeed::SetupDelta ( const std::string & databasePath)
private

◆ SetupGas()

bool Garfield::TrackHeed::SetupGas ( Medium * medium)
private

◆ SetupMaterial()

bool Garfield::TrackHeed::SetupMaterial ( Medium * medium)
private

◆ TransportDeltaElectron() [1/3]

Cluster Garfield::TrackHeed::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 )

Simulate a delta electron.

Parameters
x0,y0,z0initial position of the delta electron
t0initial time
e0initial kinetic energy of the delta electron
dx0,dy0,dz0initial direction of the delta electron

◆ TransportDeltaElectron() [2/3]

void Garfield::TrackHeed::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 )

Simulate a delta electron.

Parameters
x0,y0,z0initial position of the delta electron
t0initial time
e0initial kinetic energy of the delta electron
dx0,dy0,dz0initial direction of the delta electron
nenumber of electrons produced by the delta electron

◆ TransportDeltaElectron() [3/3]

void Garfield::TrackHeed::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 )

Simulate a delta electron.

Parameters
x0,y0,z0initial position of the delta electron
t0initial time
e0initial kinetic energy of the delta electron
dx0,dy0,dz0initial direction of the delta electron
ne,ninumber of electrons/ions produced by the delta electron

◆ TransportPhoton() [1/4]

Cluster Garfield::TrackHeed::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 )

Simulate a photon.

Parameters
x0,y0,z0initial position of the photon
t0initial time
e0initial energy of the photon
dx0,dy0,dz0initial direction of the photon

◆ TransportPhoton() [2/4]

void Garfield::TrackHeed::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 )

Simulate a photon.

Parameters
x0,y0,z0initial position of the photon
t0initial time
e0initial energy of the photon
dx0,dy0,dz0initial direction of the photon
nenumber of electrons produced by the photon

◆ TransportPhoton() [3/4]

void Garfield::TrackHeed::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 )

Simulate a photon.

Parameters
x0,y0,z0initial position of the photon
t0initial time
e0initial energy of the photon
dx0,dy0,dz0initial direction of the photon
nenumber of electrons produced by the photon
ninumber of ions produced by the photon

◆ TransportPhoton() [4/4]

void Garfield::TrackHeed::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,
int & np )

Simulate a photon.

Parameters
x0,y0,z0initial position of the photon
t0initial time
e0initial energy of the photon
dx0,dy0,dz0initial direction of the photon
nenumber of electrons produced by the photon
ninumber of ions produced by the photon
npnumber of fluorescence photons

◆ UpdateBoundingBox()

bool Garfield::TrackHeed::UpdateBoundingBox ( bool & update)
private

Member Data Documentation

◆ m_chamber

std::unique_ptr<HeedChamber> Garfield::TrackHeed::m_chamber
private

Definition at line 325 of file TrackHeed.hh.

◆ m_cluster

size_t Garfield::TrackHeed::m_cluster = 0
private

Definition at line 302 of file TrackHeed.hh.

◆ m_clusters

std::vector<Cluster> Garfield::TrackHeed::m_clusters
private

Definition at line 301 of file TrackHeed.hh.

◆ m_coulombScattering

bool Garfield::TrackHeed::m_coulombScattering = false
private

Definition at line 297 of file TrackHeed.hh.

◆ m_crossInactiveMedia

bool Garfield::TrackHeed::m_crossInactiveMedia = false
private

Definition at line 296 of file TrackHeed.hh.

◆ m_cX

double Garfield::TrackHeed::m_cX = 0.
private

Definition at line 330 of file TrackHeed.hh.

◆ m_cY

double Garfield::TrackHeed::m_cY = 0.
private

Definition at line 330 of file TrackHeed.hh.

◆ m_cZ

double Garfield::TrackHeed::m_cZ = 0.
private

Definition at line 330 of file TrackHeed.hh.

◆ m_deltaCs

std::unique_ptr<Heed::HeedDeltaElectronCS> Garfield::TrackHeed::m_deltaCs
private

Definition at line 322 of file TrackHeed.hh.

◆ m_doDeltaTransport

bool Garfield::TrackHeed::m_doDeltaTransport = true
private

Definition at line 299 of file TrackHeed.hh.

◆ m_doPhotonReabsorption

bool Garfield::TrackHeed::m_doPhotonReabsorption = true
private

Definition at line 295 of file TrackHeed.hh.

◆ m_elScat

std::unique_ptr<Heed::ElElasticScat> Garfield::TrackHeed::m_elScat
private

Definition at line 319 of file TrackHeed.hh.

◆ m_emax

double Garfield::TrackHeed::m_emax = 2.e-1
private

Definition at line 313 of file TrackHeed.hh.

◆ m_emin

double Garfield::TrackHeed::m_emin = 2.e-6
private

Definition at line 312 of file TrackHeed.hh.

◆ m_energyMesh

std::unique_ptr<Heed::EnergyMesh> Garfield::TrackHeed::m_energyMesh
private

Definition at line 315 of file TrackHeed.hh.

◆ m_fieldMap

std::unique_ptr<HeedFieldMap> Garfield::TrackHeed::m_fieldMap
private

Definition at line 326 of file TrackHeed.hh.

◆ m_gas

std::unique_ptr<Heed::GasDef> Garfield::TrackHeed::m_gas
private

Definition at line 308 of file TrackHeed.hh.

◆ m_hasActiveTrack

bool Garfield::TrackHeed::m_hasActiveTrack = false
private

Definition at line 288 of file TrackHeed.hh.

◆ m_lowSigma

std::unique_ptr<Heed::ElElasticScatLowSigma> Garfield::TrackHeed::m_lowSigma
private

Definition at line 320 of file TrackHeed.hh.

◆ m_lX

double Garfield::TrackHeed::m_lX = 0.
private

Definition at line 329 of file TrackHeed.hh.

◆ m_lY

double Garfield::TrackHeed::m_lY = 0.
private

Definition at line 329 of file TrackHeed.hh.

◆ m_lZ

double Garfield::TrackHeed::m_lZ = 0.
private

Definition at line 329 of file TrackHeed.hh.

◆ m_material

std::unique_ptr<Heed::MatterDef> Garfield::TrackHeed::m_material
private

Definition at line 309 of file TrackHeed.hh.

◆ m_matter

std::unique_ptr<Heed::HeedMatterDef> Garfield::TrackHeed::m_matter
private

Definition at line 307 of file TrackHeed.hh.

◆ m_maxStep

double Garfield::TrackHeed::m_maxStep = 100.
private

Max. step length.

Definition at line 334 of file TrackHeed.hh.

◆ m_mediumDensity

double Garfield::TrackHeed::m_mediumDensity = -1.
private

Definition at line 290 of file TrackHeed.hh.

◆ m_mediumName

std::string Garfield::TrackHeed::m_mediumName = ""
private

Definition at line 291 of file TrackHeed.hh.

◆ m_nEnergyIntervals

unsigned int Garfield::TrackHeed::m_nEnergyIntervals = 200
private

Definition at line 314 of file TrackHeed.hh.

◆ m_oneStepFly

bool Garfield::TrackHeed::m_oneStepFly = false
private

Definition at line 285 of file TrackHeed.hh.

◆ m_pairProd

std::unique_ptr<Heed::PairProd> Garfield::TrackHeed::m_pairProd
private

Definition at line 321 of file TrackHeed.hh.

◆ m_particle_def

std::unique_ptr<Heed::particle_def> Garfield::TrackHeed::m_particle_def
private

Definition at line 305 of file TrackHeed.hh.

◆ m_radStraight

double Garfield::TrackHeed::m_radStraight = 1000.
private

Bending radius beyond which to use straight-line approximation.

Definition at line 336 of file TrackHeed.hh.

◆ m_ready

bool Garfield::TrackHeed::m_ready = false
private

Definition at line 287 of file TrackHeed.hh.

◆ m_stepAngleCurved

double Garfield::TrackHeed::m_stepAngleCurved = 0.2
private

Angular step for curved lines.

Definition at line 340 of file TrackHeed.hh.

◆ m_stepAngleStraight

double Garfield::TrackHeed::m_stepAngleStraight = 0.1
private

Angular step for curved trajectories approximated by straight-line steps.

Definition at line 338 of file TrackHeed.hh.

◆ m_transferCs

std::unique_ptr<Heed::EnTransfCS> Garfield::TrackHeed::m_transferCs
private

Definition at line 318 of file TrackHeed.hh.

◆ m_useBfieldAuto

bool Garfield::TrackHeed::m_useBfieldAuto = true
private

Definition at line 298 of file TrackHeed.hh.

◆ m_usePacsOutput

bool Garfield::TrackHeed::m_usePacsOutput = false
private

Definition at line 293 of file TrackHeed.hh.


The documentation for this class was generated from the following file:
  • /builds/garfield/docs/source/Include/Garfield/TrackHeed.hh