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

Namespaces

namespace  Degrade
namespace  Magboltz
namespace  Numerics
 Collection of numerical routines.
namespace  Polygon

Classes

class  AvalancheGrid
 Calculate avalanches in a uniform electric field using avalanche statistics. More...
class  AvalancheGridSpaceCharge
 Propagates avalanches with the 2d (axi-symmetric) space-charge routine from Lippmann, Riegler (2004) in uniform background fields. More...
class  AvalancheMC
 Calculate drift lines and avalanches based on macroscopic transport coefficients, using Monte Carlo integration. More...
class  AvalancheMicroscopic
 Calculate electron drift lines and avalanches using microscopic tracking. More...
class  ComponentAnalyticField
 Semi-analytic calculation of two-dimensional configurations consisting of wires, planes, and tubes. More...
class  ComponentAnsys121
 Component for importing and interpolating two-dimensional ANSYS field maps. More...
class  ComponentAnsys123
 Component for importing and interpolating three-dimensional ANSYS field maps. More...
class  ComponentComsol
 Component for importing and interpolating Comsol field maps. More...
class  ComponentConstant
 Component with constant electric field. More...
class  ComponentCST
 Component for importing and interpolating field maps from CST. More...
class  ComponentElmer
 Component for importing field maps computed by Elmer. More...
class  ComponentElmer2d
 Component for importing two-dimensional field maps computed by Elmer. More...
class  ComponentFieldMap
 Base class for components based on finite-element field maps. More...
class  ComponentGrid
 Component for interpolating field maps on a regular mesh. More...
class  ComponentNeBem2d
 Two-dimensional implementation of the nearly exact Boundary Element Method. More...
class  ComponentNeBem3d
 Interface to neBEM. More...
class  ComponentNeBem3dMap
 Component for interpolating field maps stored in a mesh generated by neBEM. More...
class  ComponentParallelPlate
 Component for parallel-plate geometries. More...
class  ComponentTcad2d
 Interpolation in a two-dimensional field map created by Sentaurus Device. More...
class  ComponentTcad3d
 Interpolation in a three-dimensional field map created by Sentaurus Device. More...
class  ComponentTcadBase
 Interpolation in a field map created by Sentaurus Device. More...
class  ComponentUser
 Simple component with electric field given by a user function. More...
class  ComponentVoxel
 Component for interpolating field maps stored in a regular mesh. More...
class  DriftLineRKF
 Calculation of drift lines based on macroscopic transport coefficients using Runge-Kutta-Fehlberg integration. More...
class  Geometry
 Abstract base class for geometry classes. More...
class  GeometryRoot
 Use a geometry defined using the ROOT TGeo package. More...
class  GeometrySimple
 "Native" geometry, using simple shapes. More...
class  KDTree
 Main k-d tree class. More...
class  KDTreeNode
 A node in the k-d tree. More...
struct  KDTreeResult
 Search result. More...
class  MediumCdTe
 Cadmium-Telluride. More...
class  MediumConductor
 Conducting medium. More...
class  MediumDiamond
 Diamond. More...
class  MediumGaAs
 Gallium-Arsenide. More...
class  MediumGaN
 Gallium-Nitride. More...
class  MediumGas
 Base class for gas media. More...
class  MediumMagboltz
 Interface to Magboltz (version 11). More...
class  MediumPlastic
 Plastic medium. More...
class  MediumSilicon
 Solid crystalline silicon More...
class  OpticalData
 Photoabsorption cross-sections for some gases. More...
struct  Panel
 Surface panel. More...
class  PlottingEngine
 Plotting style. More...
class  QuadTree
 Quadtree search. More...
class  Random
class  RandomEngine
 Abstract base class for random number generators. More...
class  RandomEngineRoot
class  RandomEngineSTL
class  Shaper
 Class for signal processing. More...
class  Solid
 Abstract base class for solids. More...
class  SolidBox
 Box. More...
class  SolidExtrusion
 Extrusion. More...
class  SolidHole
 Box with a cylindrical hole. More...
class  SolidRidge
 Triangular prism (Toblerone bar). More...
class  SolidSphere
 Sphere. More...
class  SolidTube
 Cylindrical tube. More...
class  SolidWire
 Wire. More...
class  Track
 Abstract base class for track generation. More...
class  TrackBichsel
 Generate tracks using differential cross-sections for silicon computed by Hans Bichsel. More...
class  TrackDegrade
 Interface to Degrade. More...
class  TrackElectron
 [WIP] Ionization calculation based on MIP program (S. Biagi). More...
class  TrackHeed
 Generate tracks using Heed++. More...
class  TrackPAI
class  TrackSimple
 Generate tracks based on a cluster density given by the user. More...
class  TrackSrim
 Generate tracks based on SRIM energy loss, range and straggling tables. More...
class  TrackTrim
 Generate tracks based on TRIM output files. More...
class  Vec1Impl
class  Vec2Impl
class  Vec3Impl
class  Vector
class  ViewBase
 Base class for visualization classes. More...
class  ViewCell
 Visualize the "cell" defined in an analytic-field component. More...
class  ViewDrift
 Visualize drift lines and tracks. More...
class  ViewFEMesh
 Draw the mesh of a field-map component. More...
class  ViewField
 Visualize the potential or electric field of a component or sensor. More...
class  ViewGeometry
 Visualize a geometry defined using the "native" shapes. More...
class  ViewIsochrons
 Draw equal time contour lines. More...
class  ViewMedium
 Plot transport coefficients as function of electric and magnetic field. More...
class  ViewSignal
 Plot the signal computed by a sensor as a ROOT histogram. More...

Typedefs

typedef std::vector< std::vector< double > > KDTreeArray
using Vec3 = Vec3Impl<double>
using Vec3D = Vec3Impl<double>
using Vec3F = Vec3Impl<float>
using ViewMesh = ViewFEMesh

Enumerations

enum class  Particle {
  Electron = 0 , Ion , Hole , Positron ,
  NegativeIon , Photon
}
enum class  MPRunMode { Normal = 0 , GPUExclusive }

Functions

class GARFIELD_CLASS_NAME (Component)
 Abstract base class for components.
class GARFIELD_CLASS_NAME (Medium)
 Abstract base class for components.
void SetDefaultStyle ()
void SetSerif ()
void SetSansSerif ()
double RndmUniform ()
 Draw a random number uniformly distributed in the range [0, 1).
double RndmUniformPos ()
 Draw a random number uniformly distributed in the range (0, 1).
double RndmGaussian ()
 Draw a Gaussian random variate with mean zero and standard deviation one.
double RndmGaussian (const double mu, const double sigma)
 Draw a Gaussian random variate with mean mu and standard deviation sigma.
double RndmLorentzian (const double mu, const double gamma)
 Draw a Lorentzian random variate with mean mu and half-width at half maximum gamma.
double RndmVoigt (const double mu, const double sigma, const double gamma)
 Draw a random number according to a Voigt function with mean mu.
unsigned int RndmYuleFurry (const double mean)
 Draw a random number from a geometric distribution.
double RndmPolya (const double theta)
 Draw a Polya distributed random number.
double RndmLandau ()
 Draw a random number from a Landau distribution.
double RndmVavilov (const double rkappa, const double beta2)
 Draw a random number from a Vavilov distribution.
int RndmPoisson (const double mean)
 Draw a random number from a Poisson distribution.
double RndmHeedWF (const double w, const double f)
 Draw a random energy needed to create a single electron in a material asymptotic work function W and Fano factor F, according to Igor Smirnov's phenomenological model.
void RndmDirection (double &dx, double &dy, double &dz, const double length=1.)
 Draw a random (isotropic) direction vector.
class GARFIELD_CLASS_NAME (Sensor)
 Sensor
class GARFIELD_CLASS_NAME (TetrahedralTree)
 Helper class for searches in field maps.
void ltrim (std::string &line)
void rtrim (std::string &line)
std::vector< std::string > tokenize (const std::string &line)
bool startsWith (const std::string &line, const std::string &s)

Variables

ComponentNeBem3dgComponentNeBem3d
PlottingEngine plottingEngine

Typedef Documentation

◆ KDTreeArray

typedef std::vector<std::vector<double> > Garfield::KDTreeArray

Definition at line 19 of file KDTree.hh.

◆ Vec3

using Garfield::Vec3 = Vec3Impl<double>

Definition at line 27 of file TetrahedralTree.hh.

◆ Vec3D

using Garfield::Vec3D = Vec3Impl<double>

Definition at line 98 of file Vector.hh.

◆ Vec3F

using Garfield::Vec3F = Vec3Impl<float>

Definition at line 99 of file Vector.hh.

◆ ViewMesh

Definition at line 163 of file ViewFEMesh.hh.

Enumeration Type Documentation

◆ MPRunMode

enum class Garfield::MPRunMode
strong
Enumerator
Normal 
GPUExclusive 

Definition at line 7 of file MultiProcessInterface.hh.

◆ Particle

enum class Garfield::Particle
strong
Enumerator
Electron 
Ion 
Hole 
Positron 
NegativeIon 
Photon 

Definition at line 7 of file GarfieldConstants.hh.

Function Documentation

◆ GARFIELD_CLASS_NAME() [1/4]

class Garfield::GARFIELD_CLASS_NAME ( Component )

Abstract base class for components.

Default constructor.

Constructor

Destructor

Define the geometry.

Reset.

Get the medium at a given location (x, y, z).

Calculate the drift field at given point.

Parameters
x,y,zcoordinates [cm].
ex,ey,ezcomponents of the electric field [V/cm].
mpointer to the medium at this location.
statusstatus flag

Status flags:

        0: Inside an active medium
      > 0: Inside a wire of type X
-4 ... -1: On the side of a plane where no wires are
       -5: Inside the mesh but not in an active medium
       -6: Outside the mesh
      -10: Unknown potential type (should not occur)
    other: Other cases (should not occur)

Calculate the drift field [V/cm] and potential [V] at (x, y, z).

Calculate the drift field [V/cm] at (x, y, z).

Calculate the (drift) electrostatic potential [V] at (x, y, z).

Calculate the voltage range [V].

Calculate the weighting field at a given point and for a given electrode.

Parameters
x,y,zcoordinates [cm].
wx,wy,wzcomponents of the weighting field [1/cm].
labelname of the electrode

Calculate the weighting potential at a given point.

Parameters
x,y,zcoordinates [cm].
labelname of the electrode.
Returns
weighting potential [dimensionless].

Return the time steps at which the delayed weighting potential/field are stored/evaluated.

Calculate the delayed weighting field at a given point and time and for a given electrode.

Parameters
x,y,zcoordinates [cm].
ttime [ns].
wx,wy,wzcomponents of the weighting field [1/cm].
labelname of the electrode

Calculate the delayed weighting potential at a given point and time and for a given electrode.

Parameters
x,y,zcoordinates [cm].
ttime [ns].
labelname of the electrode

Calculate the delayed weighting potentials at a given point and for a given electrode, for a set of pre-defined times.

Calculate the magnetic field at a given point.

Parameters
x,y,zcoordinates [cm].
bx,by,bzcomponents of the magnetic field [Tesla].
statusstatus flag.

Set a constant magnetic field.

Ready for use?

Does the component have a 3D field (map)?

Get the bounding box coordinates.

Get the coordinates of the elementary cell.

Return the number of mesh elements.

Get the indices of the nodes constituting a given element.

Get the region/material of a mesh element and a flag whether it is associated to an active medium.

Return the number of mesh nodes.

Get the coordinates of a mesh node.

Integrate the normal component of the electric field over a circle.

Parameters
xc,yccentre of the circle [cm]
rradius [cm]
nInumber of intervals for the integration
Returns
charge enclosed in the circle [fC / cm]

Integrate the normal component of the electric field over a sphere.

Parameters
xc,yc,zccentre of the sphere [cm]
rradius of the sphere [cm]
nInumber of integration intervals in phi and theta
Returns
charge enclosed in the sphere [fC]

Integrate the normal component of the electric field over a parallelogram.

Parameters
x0,y0,z0coordinates of one of the corners [cm]
dx1,dy1,dz1vector to one of the adjacent corners [cm]
dx2,dy2,dz2vector to the other adjacent corner [cm]
nU,nVnumber of integration points in the two directions
Returns
flux [V cm]

Integrate the normal component of the weighting field over a parallelogram.

Integrate the electric field flux through a line from (x0,y0,z0) to (x1,y1,z1) along a direction (xp,yp,zp).

Parameters
x0,y0,z0coordinates of the starting point
x1,y1,z1coordinates of the end point
xp,yp,zpnormal vector
nInumber of intervals for the integration
isigninclude both negative and positive contributions (0) or only contributions with a given polarity (+1,-1)

Determine whether the line between two points crosses a wire.

Parameters
x0,y0,z0first point [cm].
x1,y1,z1second point [cm]
xc,yc,zcpoint [cm] where the line crosses the wire or the coordinates of the wire centre.
centreflag whether to return the coordinates of the line-wire crossing point or of the wire centre.
rcradius [cm] of the wire.

Determine whether a particle is inside the trap radius of a wire.

Parameters
q0charge of the particle [in elementary charges].
x0,y0,z0position [cm] of the particle.
xw,ywcoordinates of the wire (if applicable).
rwradius of the wire (if applicable).

Determine whether the line between two points crosses a plane.

Enable simple periodicity in the $x$ direction.

Enable simple periodicity in the $y$ direction.

Enable simple periodicity in the $z$ direction.

Return periodicity flags.

Enable mirror periodicity in the $x$ direction.

Enable mirror periodicity in the $y$ direction.

Enable mirror periodicity in the $y$ direction.

Return mirror periodicity flags.

Enable axial periodicity in the $x$ direction.

Enable axial periodicity in the $y$ direction.

Enable axial periodicity in the $z$ direction.

Return axial periodicity flags.

Enable rotation symmetry around the $x$ axis.

Enable rotation symmetry around the $y$ axis.

Enable rotation symmetry around the $z$ axis.

Return rotation symmetry flags.

Enable triangular periodicity in the $xy$ plane.

Enable triangular periodicity in the $xz$ plane.

Enable triangular periodicity in the $yz$ plane.

Switch on debugging messages.

Switch off debugging messages.

Does the component have a non-zero magnetic field?

Does the component have maps of the Townsend coefficient?

Does the component have maps of the attachment coefficient?

Does the component have maps of the low-field mobility?

Does the component have velocity maps?

Get the electron attachment coefficient.

Get the hole attachment coefficient.

Get the hole Mobility coefficient.

Get the electron Townsend coefficient.

Get the hole Townsend coefficient.

Get the electron drift velocity.

Get the hole drift velocity.

Create and initialise GPU Transfer class

Class name.

Pointer to the geometry.

Constant magnetic field.

Ready for use?

Switch on/off debugging messages

Simple periodicity in x, y, z.

Mirror periodicity in x, y, z.

Axial periodicity in x, y, z.

Rotation symmetry around x-axis, y-axis, z-axis.

Triangle symmetric in the xy, xz, and yz plane.

Triangle symmetric octant of imported map (0 < phi < Pi/4 --> octant 1).

Octants where |x| >= |y|

Time steps at which the delayed weighting potentials/fields are stored.

Reset the component.

Verify periodicities.

Definition at line 1 of file Component.hh.

31 {
32 public:
33#ifdef __GPUCOMPILE__
34 GARFIELD_CLASS_NAME(Component)() = default;
35#else
37 GARFIELD_CLASS_NAME(Component)() = delete;
39 GARFIELD_CLASS_NAME(Component)(const std::string& name);
40#endif
42 virtual ~GARFIELD_CLASS_NAME(Component)() {};
43
44#ifndef __GPUCOMPILE__
46 virtual void SetGeometry(Geometry* geo);
48 virtual void Clear();
49
51 virtual Medium* GetMedium(const double x, const double y, const double z);
52
70 virtual void ElectricField(const double x, const double y, const double z,
71 double& ex, double& ey, double& ez, Medium*& m,
72 int& status) = 0;
74 virtual void ElectricField(const double x, const double y, const double z,
75 double& ex, double& ey, double& ez, double& v,
76 Medium*& m, int& status) = 0;
77#else
78 __device__ void ElectricField(const cuda_t xin, const cuda_t yin,
79 const cuda_t zin, cuda_t& ex, cuda_t& ey,
80 cuda_t& ez, MediumGPU*& m, int& status);
81#endif
82
83#ifndef __GPUCOMPILE__
85 std::array<double, 3> ElectricField(const double x, const double y,
86 const double z);
88 virtual double ElectricPotential(const double x, const double y,
89 const double z);
91 virtual bool GetVoltageRange(double& vmin, double& vmax) = 0;
92
98 virtual void WeightingField(const double x, const double y, const double z,
99 double& wx, double& wy, double& wz,
100 const std::string& label);
106 virtual double WeightingPotential(const double x, const double y,
107 const double z, const std::string& label);
108
112 virtual const std::vector<double>& DelayedSignalTimes(
113 const std::string& /*label*/) {
114 return m_wdtimes;
115 }
123 virtual void DelayedWeightingField(const double x, const double y,
124 const double z, const double t, double& wx,
125 double& wy, double& wz,
126 const std::string& label);
127
134 virtual double DelayedWeightingPotential(const double x, const double y,
135 const double z, const double t,
136 const std::string& label);
137
141 virtual void DelayedWeightingPotentials(const double x, const double y,
142 const double z,
143 const std::string& label,
144 std::vector<double>& dwp);
151 virtual void MagneticField(const double x, const double y, const double z,
152 double& bx, double& by, double& bz, int& status);
154 void SetMagneticField(const double bx, const double by, const double bz);
155
157 virtual bool IsReady() { return m_ready; }
159 virtual bool Is3d() { return true; }
160
162 virtual bool GetBoundingBox(double& xmin, double& ymin, double& zmin,
163 double& xmax, double& ymax, double& zmax);
164
166 virtual bool GetElementaryCell(double& xmin, double& ymin, double& zmin,
167 double& xmax, double& ymax, double& zmax);
168 // Get the x-length of the elementary cell.
169 double CellSizeX();
170 // Get the y-length of the elementary cell.
171 double CellSizeY();
172 // Get the z-length of the elementary cell.
173 double CellSizeZ();
174
176 virtual size_t GetNumberOfElements() const { return 0; }
178 virtual bool GetElementNodes(const size_t /*i*/,
179 std::vector<size_t>& /*nodes*/) const {
180 return false;
181 }
184 virtual bool GetElementRegion(const size_t /*i*/, size_t& /*mat*/,
185 bool& /*drift*/) const {
186 return false;
187 }
189 virtual size_t GetNumberOfNodes() const { return 0; }
191 virtual bool GetNode(const size_t i, double& x, double& y, double& z) const;
192
200 double IntegrateFluxCircle(const double xc, const double yc, const double r,
201 const unsigned int nI = 50);
209 double IntegrateFluxSphere(const double xc, const double yc, const double zc,
210 const double r, const unsigned int nI = 20);
211
220 double IntegrateFluxParallelogram(
221 const double x0, const double y0, const double z0, const double dx1,
222 const double dy1, const double dz1, const double dx2, const double dy2,
223 const double dz2, const unsigned int nU = 20, const unsigned int nV = 20);
224
227 double IntegrateWeightingFluxParallelogram(
228 const std::string& label, const double x0, const double y0,
229 const double z0, const double dx1, const double dy1, const double dz1,
230 const double dx2, const double dy2, const double dz2,
231 const unsigned int nU = 20, const unsigned int nV = 20);
232
242 double IntegrateFluxLine(const double x0, const double y0, const double z0,
243 const double x1, const double y1, const double z1,
244 const double xp, const double yp, const double zp,
245 const unsigned int nI, const int isign = 0);
246
256 virtual bool CrossedWire(const double x0, const double y0, const double z0,
257 const double x1, const double y1, const double z1,
258 double& xc, double& yc, double& zc,
259 const bool centre, double& rc);
266 virtual bool InTrapRadius(const double q0, const double x0, const double y0,
267 const double z0, double& xw, double& yw,
268 double& rw);
271 virtual bool CrossedPlane(const double x0, const double y0, const double z0,
272 const double x1, const double y1, const double z1,
273 double& xc, double& yc, double& zc);
274
276 void EnablePeriodicityX(const bool on = true) {
277 m_periodic[0] = on;
278 UpdatePeriodicity();
279 }
281 void EnablePeriodicityY(const bool on = true) {
282 m_periodic[1] = on;
283 UpdatePeriodicity();
284 }
286 void EnablePeriodicityZ(const bool on = true) {
287 m_periodic[2] = on;
288 UpdatePeriodicity();
289 }
291 void IsPeriodic(bool& perx, bool& pery, bool& perz) {
292 perx = m_periodic[0];
293 pery = m_periodic[1];
294 perz = m_periodic[2];
295 }
296
298 void EnableMirrorPeriodicityX(const bool on = true) {
299 m_mirrorPeriodic[0] = on;
300 UpdatePeriodicity();
301 }
303 void EnableMirrorPeriodicityY(const bool on = true) {
304 m_mirrorPeriodic[1] = on;
305 UpdatePeriodicity();
306 }
308 void EnableMirrorPeriodicityZ(const bool on = true) {
309 m_mirrorPeriodic[2] = on;
310 UpdatePeriodicity();
311 }
313 void IsMirrorPeriodic(bool& perx, bool& pery, bool& perz) {
314 perx = m_mirrorPeriodic[0];
315 pery = m_mirrorPeriodic[1];
316 perz = m_mirrorPeriodic[2];
317 }
318
320 void EnableAxialPeriodicityX(const bool on = true) {
321 m_axiallyPeriodic[0] = on;
322 UpdatePeriodicity();
323 }
325 void EnableAxialPeriodicityY(const bool on = true) {
326 m_axiallyPeriodic[1] = on;
327 UpdatePeriodicity();
328 }
330 void EnableAxialPeriodicityZ(const bool on = true) {
331 m_axiallyPeriodic[2] = on;
332 UpdatePeriodicity();
333 }
335 void IsAxiallyPeriodic(bool& perx, bool& pery, bool& perz) {
336 perx = m_axiallyPeriodic[0];
337 pery = m_axiallyPeriodic[1];
338 perz = m_axiallyPeriodic[2];
339 }
340
342 void EnableRotationSymmetryX(const bool on = true) {
343 m_rotationSymmetric[0] = on;
344 UpdatePeriodicity();
345 }
347 void EnableRotationSymmetryY(const bool on = true) {
348 m_rotationSymmetric[1] = on;
349 UpdatePeriodicity();
350 }
352 void EnableRotationSymmetryZ(const bool on = true) {
353 m_rotationSymmetric[2] = on;
354 UpdatePeriodicity();
355 }
357 void IsRotationSymmetric(bool& rotx, bool& roty, bool& rotz) {
358 rotx = m_rotationSymmetric[0];
359 roty = m_rotationSymmetric[1];
360 rotz = m_rotationSymmetric[2];
361 }
362
364 void EnableTriangleSymmetricXY(const bool on = true,
365 const bool oct = 2) {
366 m_triangleSymmetric[0] = on;
367 m_triangleSymmetricOct = oct;
368 m_mirrorPeriodic[0] = on;
369 m_mirrorPeriodic[1] = on;
370 }
372 void EnableTriangleSymmetricXZ(const bool on = true,
373 const bool oct = 2) {
374 m_triangleSymmetric[1] = on;
375 m_triangleSymmetricOct = oct;
376 m_mirrorPeriodic[0] = on;
377 m_mirrorPeriodic[2] = on;
378 }
380 void EnableTriangleSymmetricYZ(const bool on = true,
381 const bool oct = 2) {
382 m_triangleSymmetric[2] = on;
383 m_triangleSymmetricOct = oct;
384 m_mirrorPeriodic[1] = on;
385 m_mirrorPeriodic[2] = on;
386 }
387
389 void EnableDebugging(const bool on = true) { m_debug = on; }
391 void DisableDebugging() { m_debug = false; }
392
394 virtual bool HasMagneticField() const;
395
397 virtual bool HasTownsendMap() const { return false; }
399 virtual bool HasAttachmentMap() const { return false; }
401 virtual bool HasMobilityMap() const { return false; }
403 virtual bool HasVelocityMap() const { return false; }
404
406 virtual bool ElectronAttachment(const double /*x*/, const double /*y*/,
407 const double /*z*/, double& eta) {
408 eta = 0;
409 return false;
410 }
412 virtual bool HoleAttachment(const double /*x*/, const double /*y*/,
413 const double /*z*/, double& eta) {
414 eta = 0;
415 return false;
416 }
417
418 // Get the electron mobility coefficient.
419 virtual bool ElectronMobility(const double /*x*/, const double /*y*/,
420 const double /*z*/, double& mu) {
421 mu = 0;
422 return false;
423 }
425 virtual bool HoleMobility(const double /*x*/, const double /*y*/,
426 const double /*z*/, double& mu) {
427 mu = 0;
428 return false;
429 }
430
432 virtual bool ElectronTownsend(const double /*x*/, const double /*y*/,
433 const double /*z*/, double& alpha) {
434 alpha = 0;
435 return false;
436 }
438 virtual bool HoleTownsend(const double /*x*/, const double /*y*/,
439 const double /*z*/, double& alpha) {
440 alpha = 0;
441 return false;
442 }
444 virtual bool ElectronVelocity(const double /*x*/, const double /*y*/,
445 const double /*z*/, double& vx, double& vy,
446 double& vz) {
447 vx = vy = vz = 0;
448 return false;
449 }
451 virtual bool HoleVelocity(const double /*x*/, const double /*y*/,
452 const double /*z*/, double& vx, double& vy,
453 double& vz) {
454 vx = vy = vz = 0;
455 return false;
456 }
457
458 virtual double StepSizeHint() { return -1.; }
459
461 virtual double CreateGPUTransferObject(ComponentGPU*& comp_gpu);
462
463 protected:
465 std::string m_className = "Component";
466
468 Geometry* m_geometry = nullptr;
469
471 std::array<double, 3> m_b0 = {{0., 0., 0.}};
472#endif
474 bool m_ready = false;
475
476#ifdef __GPUCOMPILE__
478 bool m_periodic[3] = {false, false, false};
480 bool m_mirrorPeriodic[3] = {false, false, false};
482 bool m_axiallyPeriodic[3] = {false, false, false};
484 bool m_rotationSymmetric[3] = {false, false, false};
486 bool m_triangleSymmetric[3] = {false, false, false};
488 int m_triangleSymmetricOct = 0;
490 const int m_triangleOctRules[4] = {1, 4, 5, 8};
491 bool m_outsideCone = false;
492#else
494 bool m_debug = false;
495
497 std::array<bool, 3> m_periodic = {{false, false, false}};
499 std::array<bool, 3> m_mirrorPeriodic = {{false, false, false}};
501 std::array<bool, 3> m_axiallyPeriodic = {{false, false, false}};
503 std::array<bool, 3> m_rotationSymmetric = {{false, false, false}};
505 std::array<bool, 3> m_triangleSymmetric = {{false, false, false}};
507 int m_triangleSymmetricOct = 0;
509 const std::array<int, 4> m_triangleOctRules = {1, 4, 5, 8};
510 bool m_outsideCone = false;
511#endif
512
513#ifndef __GPUCOMPILE__
515 std::vector<double> m_wdtimes;
516
518 virtual void Reset() = 0;
520 virtual void UpdatePeriodicity() = 0;
521
522 private:
523 double IntegrateFluxParallelogram(const double x0, const double y0,
524 const double z0, const double dx1,
525 const double dy1, const double dz1,
526 const double dx2, const double dy2,
527 const double dz2, const unsigned int nU,
528 const unsigned int nV, const bool wfield,
529 const std::string& label);
530#else
531
532// include parts from derived class due to big performance hit from using
533// virtual methods
534// TODO GPU TN: It isn't clear to me that we actually need (at least) the Ansys
535// include - all the code is protected by an ifndef __GPUCOMPILE__, whereas it
536// will be defined when these includes are made
537#include "ComponentAnsys123.hh"
538#include "ComponentFieldMap.hh"
539
540 friend class ComponentAnsys123;
541 friend class ComponentComsol;
542 friend class ComponentElmer;
543 friend class ComponentFieldMap;
544 friend class Component;
545
546 // enum to mimic polymorphism
547 enum class ComponentType {
548 Component = 0,
553 };
554
555 ComponentType m_ComponentType{ComponentType::Component};
556
557#endif
558};
#define GARFIELD_CLASS_NAME(name)
#define __device__
Definition Vector.hh:12
Component for importing and interpolating three-dimensional ANSYS field maps.
Component for importing and interpolating Comsol field maps.
Component for importing field maps computed by Elmer.
Base class for components based on finite-element field maps.
Abstract base class for geometry classes.
Definition Geometry.hh:13
class GARFIELD_CLASS_NAME(Component)
Abstract base class for components.
Definition Component.hh:31

◆ GARFIELD_CLASS_NAME() [2/4]

class Garfield::GARFIELD_CLASS_NAME ( Medium )

Abstract base class for components.

Constructor

Destructor

Return the id number of the class instance.

Get the medium name/identifier.

Is this medium a gas?

Is this medium a semiconductor?

Is this medium a conductor?

Set the temperature [K].

Get the temperature [K].

Set the relative static dielectric constant.

Get the relative static dielectric constant.

Get number of components of the medium.

Get the name and fraction of a given component.

Set the effective atomic number.

Get the effective atomic number.

Set the effective atomic weight.

Get the effective atomic weight.

Set the number density [cm-3].

Get the number density [cm-3].

Set the mass density [g/cm3].

Get the mass density [g/cm3].

Switch electron/ion/hole transport on/off.

Make the medium ionisable or non-ionisable.

Is charge carrier transport enabled in this medium?

Does the medium have electron scattering rates?

Is charge deposition by charged particles/photon enabled in this medium?

Set the W value (average energy to produce an electron/ion or e/h pair).

Get the W value.

Set the Fano factor.

Get the Fano factor.

Plot the drift velocity as function of the electric field.

Plot the diffusion coefficients as function of the electric field.

Plot the Townsend coefficient(s) as function of the electric field.

Plot the attachment coefficient(s) as function of the electric field.

Plot Townsend and attachment coefficients.

Drift velocity [cm / ns]

Flux (mean velocity; shorthand: wv) and bulk (center of mass velocity; shorthand: wr) drift velocity [cm / ns]

Longitudinal and transverse diffusion coefficients [cm1/2]

Diffusion tensor: diagonal elements are the diffusion coefficients [cm] along e, btrans, e x b, off-diagonal elements are the covariances

Ionisation coefficient [cm-1]

Attachment coefficient [cm-1]

TOF Ionisation rate [ns-1]

TOF Attachment Rate [ns-1]

Lorentz angle

Low-field mobility [cm2 V-1 ns-1]

Dispersion relation (energy vs. wave vector)

Sample the momentum vector for a given energy (only meaningful in semiconductors).

Null-collision rate [ns-1]

Collision rate [ns-1] for given electron energy

Sample the collision type. Update energy and direction vector.

Drift velocity [cm / ns]

Longitudinal and transverse diffusion coefficients [cm1/2]

Diffusion tensor

Ionisation coefficient [cm-1]

Attachment coefficient [cm-1]

Low-field mobility [cm2 V-1 ns-1]

Ion drift velocity [cm / ns]

Longitudinal and transverse diffusion coefficients [cm1/2]

Dissociation coefficient

Low-field ion mobility [cm2 V-1 ns-1]

Negative ion drift velocity [cm / ns]

Low-field negative ion mobility [cm2 V-1 ns-1]

Set the range of fields to be covered by the transport tables.

Set the fields and E-B angles to be used in the transport tables.

Get the fields and E-B angles used in the transport tables.

Set an entry in the table of drift speeds along E.

Get an entry in the table of drift speeds along E.

Set an entry in the table of drift speeds along ExB.

Get an entry in the table of drift speeds along ExB.

Set an entry in the table of drift speeds along Btrans.

Get an entry in the table of drift speeds along Btrans.

Set an entry in the table of flux drift speeds.

Get an entry in the table of flux drift speeds.

Set an entry in the table of bulk drift speeds.

Get an entry in the table of bulk drift speeds.

Set an entry in the table of longitudinal diffusion coefficients.

Get an entry in the table of longitudinal diffusion coefficients.

Set an entry in the table of transverse diffusion coefficients.

Get an entry in the table of transverse diffusion coefficients.

Set an entry in the table of Townsend coefficients.

Get an entry in the table of Townsend coefficients.

Set an entry in the table of attachment coefficients.

Get an entry in the table of attachment coefficients.

Set an entry in the table of ionization rate of TOF.

Get an entry in the table of ionization rate of TOF.

Set an entry in the table of attachment rate of TOF.

Get an entry in the table of attachment rate of TOF.

Set an entry in the table of Lorentz angles.

Get an entry in the table of Lorentz angles.

Set an entry in the table of drift speeds along E.

Get an entry in the table of drift speeds along E.

Set an entry in the table of drift speeds along ExB.

Get an entry in the table of drift speeds along ExB.

Set an entry in the table of drift speeds along Btrans.

Get an entry in the table of drift speeds along Btrans.

Set an entry in the table of longitudinal diffusion coefficients.

Get an entry in the table of longitudinal diffusion coefficients.

Set an entry in the table of transverse diffusion coefficients.

Get an entry in the table of transverse diffusion coefficients.

Set an entry in the table of Townsend coefficients.

Get an entry in the table of Townsend coefficients.

Set an entry in the table of attachment coefficients.

Get an entry in the table of attachment coefficients.

Initialise the table of ion mobilities from a list of electric fields and corresponding mobilities. The mobilities will be interpolated at the electric fields of the currently set grid.

Set an entry in the table of ion mobilities.

Get an entry in the table of ion mobilities.

Set an entry in the table of longitudinal diffusion coefficients.

Get an entry in the table of longitudinal diffusion coefficients.

Set an entry in the table of transverse diffusion coefficients.

Get an entry in the table of transverse diffusion coefficients.

Set an entry in the table of dissociation coefficients.

Get an entry in the table of dissociation coefficients.

Set an entry in the table of negative ion mobilities.

Get an entry in the table of negative ion mobilities.

Reset all tables of transport parameters.

Select the extrapolation method for fields below/above the table range. Possible options are "constant", "linear", and "exponential".

Set the degree of polynomial interpolation (usually 2).

Get the energy range [eV] of the available optical data.

Get the complex dielectric function at a given energy.

Switch on/off debugging messages

Create and initialise GPU Transfer class

Definition at line 1 of file Medium.hh.

35 {
36 public:
37#ifdef __GPUCOMPILE__
39 GARFIELD_CLASS_NAME(Medium)() = default;
41 ~GARFIELD_CLASS_NAME(Medium)() {};
42#else
44 GARFIELD_CLASS_NAME(Medium)();
46 virtual ~GARFIELD_CLASS_NAME(Medium)();
47#endif
48
50 __DEVICE__ int GetId() const { return m_id; }
51
52#ifndef __GPUCOMPILE__
54 const std::string& GetName() const { return m_name; }
56 virtual bool IsGas() const { return false; }
58 virtual bool IsSemiconductor() const { return false; }
60 virtual bool IsConductor() const { return false; }
61
63 void SetTemperature(const double t);
65 double GetTemperature() const { return m_temperature; }
66 // Set the pressure [Torr].
67 void SetPressure(const double p);
68 // Get the pressure [Torr].
69 double GetPressure() const { return m_pressure; }
71 void SetDielectricConstant(const double eps);
73 double GetDielectricConstant() const { return m_epsilon; }
74
76 unsigned int GetNumberOfComponents() const { return m_nComponents; }
78 virtual void GetComponent(const unsigned int i, std::string& label,
79 double& f);
81 virtual void SetAtomicNumber(const double z);
83 virtual double GetAtomicNumber() const { return m_z; }
85 virtual void SetAtomicWeight(const double a);
87 virtual double GetAtomicWeight() const { return m_a; }
89 virtual void SetNumberDensity(const double n);
91 virtual double GetNumberDensity() const { return m_density; }
93 virtual void SetMassDensity(const double rho);
95 virtual double GetMassDensity() const;
96
98 virtual void EnableDrift(const bool on = true) { m_driftable = on; }
100 virtual void EnablePrimaryIonisation(const bool on = true) {
101 m_ionisable = on;
102 }
103#endif
104
106 __DEVICE__ bool IsDriftable() const { return m_driftable; }
108 __DEVICE__ bool IsMicroscopic() const { return m_microscopic; }
109
110#ifndef __GPUCOMPILE__
112 bool IsIonisable() const { return m_ionisable; }
113
115 void SetW(const double w) { m_w = w; }
117 double GetW() const { return m_w; }
119 void SetFanoFactor(const double f) { m_fano = f; }
121 double GetFanoFactor() const { return m_fano; }
122
124 void PlotVelocity(const std::string& carriers, TPad* pad);
126 void PlotDiffusion(const std::string& carriers, TPad* pad);
128 void PlotTownsend(const std::string& carriers, TPad* pad);
130 void PlotAttachment(const std::string& carriers, TPad* pad);
132 void PlotAlphaEta(const std::string& carriers, TPad* pad);
133
134 // Transport parameters for electrons
136 virtual bool ElectronVelocity(const double ex, const double ey,
137 const double ez, const double bx,
138 const double by, const double bz, double& vx,
139 double& vy, double& vz);
142 virtual bool ElectronVelocityFluxBulk(const double ex, const double ey,
143 const double ez, const double bx,
144 const double by, const double bz,
145 double& wv, double& wr);
147 virtual bool ElectronDiffusion(const double ex, const double ey,
148 const double ez, const double bx,
149 const double by, const double bz, double& dl,
150 double& dt);
154 virtual bool ElectronDiffusion(const double ex, const double ey,
155 const double ez, const double bx,
156 const double by, const double bz,
157 double cov[3][3]);
159 virtual bool ElectronTownsend(const double ex, const double ey,
160 const double ez, const double bx,
161 const double by, const double bz,
162 double& alpha);
164 virtual bool ElectronAttachment(const double ex, const double ey,
165 const double ez, const double bx,
166 const double by, const double bz,
167 double& eta);
169 virtual bool ElectronTOFIonisation(const double ex, const double ey,
170 const double ez, const double bx,
171 const double by, const double bz,
172 double& riontof);
174 virtual bool ElectronTOFAttachment(const double ex, const double ey,
175 const double ez, const double bx,
176 const double by, const double bz,
177 double& ratttof);
179 virtual bool ElectronLorentzAngle(const double ex, const double ey,
180 const double ez, const double bx,
181 const double by, const double bz,
182 double& lor);
184 virtual double ElectronMobility();
185 // Microscopic electron transport properties
186
188 virtual double GetElectronEnergy(const double px, const double py,
189 const double pz, double& vx, double& vy,
190 double& vz, const int band = 0);
193 virtual void GetElectronMomentum(const double e, double& px, double& py,
194 double& pz, int& band);
195
197 virtual double GetElectronNullCollisionRate(const int band = 0);
198#endif
199
200#ifdef __GPUCOMPILE__
201
202 __device__ cuda_t GetElectronCollisionRate(const cuda_t e, const int band);
203
204 __device__ bool ElectronCollision(const cuda_t e, int& type, int& level,
205 cuda_t& e1, cuda_t& dx, cuda_t& dy,
206 cuda_t& dz, Particle* secondaries_type,
207 cuda_t* secondaries_energy,
208 int& num_secondaries, int& ndxc, int& band);
209#else
211 virtual double GetElectronCollisionRate(const double e, const int band = 0);
212 struct Secondary {
213 Particle type = Particle::Electron;
214 double energy = 0.;
215 double time = 0.;
216 double distance = 0.;
217 };
219 virtual bool ElectronCollision(const double e, int& type, int& level,
220 double& e1, double& dx, double& dy, double& dz,
221 std::vector<Secondary>& secondaries,
222 int& band);
223#endif
224
225#ifndef __GPUCOMPILE__
226 // Transport parameters for holes
228 virtual bool HoleVelocity(const double ex, const double ey, const double ez,
229 const double bx, const double by, const double bz,
230 double& vx, double& vy, double& vz);
232 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
233 const double bx, const double by, const double bz,
234 double& dl, double& dt);
236 virtual bool HoleDiffusion(const double ex, const double ey, const double ez,
237 const double bx, const double by, const double bz,
238 double cov[3][3]);
240 virtual bool HoleTownsend(const double ex, const double ey, const double ez,
241 const double bx, const double by, const double bz,
242 double& alpha);
244 virtual bool HoleAttachment(const double ex, const double ey, const double ez,
245 const double bx, const double by, const double bz,
246 double& eta);
248 virtual double HoleMobility();
249
250 // Transport parameters for ions
252 virtual bool IonVelocity(const double ex, const double ey, const double ez,
253 const double bx, const double by, const double bz,
254 double& vx, double& vy, double& vz);
255 bool HasIonVelocity() const { return !(m_iVel.empty() && m_iMob.empty()); }
257 virtual bool IonDiffusion(const double ex, const double ey, const double ez,
258 const double bx, const double by, const double bz,
259 double& dl, double& dt);
261 virtual bool IonDissociation(const double ex, const double ey,
262 const double ez, const double bx,
263 const double by, const double bz, double& diss);
265 virtual double IonMobility();
266
268 virtual bool NegativeIonVelocity(const double ex, const double ey,
269 const double ez, const double bx,
270 const double by, const double bz, double& vx,
271 double& vy, double& vz);
273 virtual double NegativeIonMobility();
274
276 void SetFieldGrid(double emin, double emax, const size_t ne, bool logE,
277 double bmin = 0., double bmax = 0., const size_t nb = 1,
278 double amin = HalfPi, double amax = HalfPi,
279 const size_t na = 1);
281 void SetFieldGrid(const std::vector<double>& efields,
282 const std::vector<double>& bfields,
283 const std::vector<double>& angles);
285 void GetFieldGrid(std::vector<double>& efields, std::vector<double>& bfields,
286 std::vector<double>& angles);
287
289 bool SetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia,
290 const double v) {
291 return SetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
292 }
294 bool GetElectronVelocityE(const size_t ie, const size_t ib, const size_t ia,
295 double& v) {
296 return GetEntry(ie, ib, ia, "ElectronVelocityE", m_eVelE, v);
297 }
299 bool SetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia,
300 const double v) {
301 return SetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
302 }
304 bool GetElectronVelocityExB(const size_t ie, const size_t ib, const size_t ia,
305 double& v) {
306 return GetEntry(ie, ib, ia, "ElectronVelocityExB", m_eVelX, v);
307 }
309 bool SetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia,
310 const double v) {
311 return SetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
312 }
314 bool GetElectronVelocityB(const size_t ie, const size_t ib, const size_t ia,
315 double& v) {
316 return GetEntry(ie, ib, ia, "ElectronVelocityB", m_eVelB, v);
317 }
319 bool SetElectronFluxVelocity(const size_t ie, const size_t ib,
320 const size_t ia, const double v) {
321 return SetEntry(ie, ib, ia, "ElectronFluxVelocity", m_eVelWv, v);
322 }
324 bool GetElectronFluxVelocity(const size_t ie, const size_t ib,
325 const size_t ia, double& v) {
326 return GetEntry(ie, ib, ia, "ElectronFluxVelocity", m_eVelWv, v);
327 }
329 bool SetElectronBulkVelocity(const size_t ie, const size_t ib,
330 const size_t ia, const double v) {
331 return SetEntry(ie, ib, ia, "ElectronBulkVelocity", m_eVelWr, v);
332 }
334 bool GetElectronBulkVelocity(const size_t ie, const size_t ib,
335 const size_t ia, double& v) {
336 return GetEntry(ie, ib, ia, "ElectronBulkVelocity", m_eVelWr, v);
337 }
339 bool SetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
340 const size_t ia, const double dl) {
341 return SetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
342 }
344 bool GetElectronLongitudinalDiffusion(const size_t ie, const size_t ib,
345 const size_t ia, double& dl) {
346 return GetEntry(ie, ib, ia, "ElectronLongitudinalDiffusion", m_eDifL, dl);
347 }
349 bool SetElectronTransverseDiffusion(const size_t ie, const size_t ib,
350 const size_t ia, const double dt) {
351 return SetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
352 }
354 bool GetElectronTransverseDiffusion(const size_t ie, const size_t ib,
355 const size_t ia, double& dt) {
356 return GetEntry(ie, ib, ia, "ElectronTransverseDiffusion", m_eDifT, dt);
357 }
359 bool SetElectronTownsend(const size_t ie, const size_t ib, const size_t ia,
360 const double alpha) {
361 return SetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
362 }
364 bool GetElectronTownsend(const size_t ie, const size_t ib, const size_t ia,
365 double& alpha) {
366 return GetEntry(ie, ib, ia, "ElectronTownsend", m_eAlp, alpha);
367 }
369 bool SetElectronAttachment(const size_t ie, const size_t ib, const size_t ia,
370 const double eta) {
371 return SetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
372 }
374 bool GetElectronAttachment(const size_t ie, const size_t ib, const size_t ia,
375 double& eta) {
376 return GetEntry(ie, ib, ia, "ElectronAttachment", m_eAtt, eta);
377 }
379 bool SetElectronTOFIonisation(const size_t ie, const size_t ib,
380 const size_t ia, const double v) {
381 return SetEntry(ie, ib, ia, "ElectronTOFIonisation", m_eRIon, v);
382 }
384 bool GetElectronTOFIonisation(const size_t ie, const size_t ib,
385 const size_t ia, double& v) {
386 return GetEntry(ie, ib, ia, "ElectronTOFIonisation", m_eRIon, v);
387 }
389 bool SetElectronTOFAttachment(const size_t ie, const size_t ib,
390 const size_t ia, const double v) {
391 return SetEntry(ie, ib, ia, "ElectronTOFAttachment", m_eRAtt, v);
392 }
394 bool GetElectronTOFAttachment(const size_t ie, const size_t ib,
395 const size_t ia, double& v) {
396 return GetEntry(ie, ib, ia, "ElectronTOFAttachment", m_eRAtt, v);
397 }
398
400 bool SetElectronLorentzAngle(const size_t ie, const size_t ib,
401 const size_t ia, const double lor) {
402 return SetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
403 }
405 bool GetElectronLorentzAngle(const size_t ie, const size_t ib,
406 const size_t ia, double& lor) {
407 return GetEntry(ie, ib, ia, "ElectronLorentzAngle", m_eLor, lor);
408 }
409
411 bool SetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia,
412 const double v) {
413 return SetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
414 }
416 bool GetHoleVelocityE(const size_t ie, const size_t ib, const size_t ia,
417 double& v) {
418 return GetEntry(ie, ib, ia, "HoleVelocityE", m_hVelE, v);
419 }
421 bool SetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia,
422 const double v) {
423 return SetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
424 }
426 bool GetHoleVelocityExB(const size_t ie, const size_t ib, const size_t ia,
427 double& v) {
428 return GetEntry(ie, ib, ia, "HoleVelocityExB", m_hVelX, v);
429 }
431 bool SetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia,
432 const double v) {
433 return SetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
434 }
436 bool GetHoleVelocityB(const size_t ie, const size_t ib, const size_t ia,
437 double& v) {
438 return GetEntry(ie, ib, ia, "HoleVelocityB", m_hVelB, v);
439 }
440
442 bool SetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
443 const size_t ia, const double dl) {
444 return SetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
445 }
447 bool GetHoleLongitudinalDiffusion(const size_t ie, const size_t ib,
448 const size_t ia, double& dl) {
449 return GetEntry(ie, ib, ia, "HoleLongitudinalDiffusion", m_hDifL, dl);
450 }
452 bool SetHoleTransverseDiffusion(const size_t ie, const size_t ib,
453 const size_t ia, const double dt) {
454 return SetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
455 }
457 bool GetHoleTransverseDiffusion(const size_t ie, const size_t ib,
458 const size_t ia, double& dt) {
459 return GetEntry(ie, ib, ia, "HoleTransverseDiffusion", m_hDifT, dt);
460 }
462 bool SetHoleTownsend(const size_t ie, const size_t ib, const size_t ia,
463 const double alpha) {
464 return SetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
465 }
467 bool GetHoleTownsend(const size_t ie, const size_t ib, const size_t ia,
468 double& alpha) {
469 return GetEntry(ie, ib, ia, "HoleTownsend", m_hAlp, alpha);
470 }
472 bool SetHoleAttachment(const size_t ie, const size_t ib, const size_t ia,
473 const double eta) {
474 return SetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
475 }
477 bool GetHoleAttachment(const size_t ie, const size_t ib, const size_t ia,
478 double& eta) {
479 return GetEntry(ie, ib, ia, "HoleAttachment", m_hAtt, eta);
480 }
481
486 bool SetIonMobility(const std::vector<double>& fields,
487 const std::vector<double>& mobilities,
488 const bool negativeIons = false);
490 bool SetIonMobility(const size_t ie, const size_t ib, const size_t ia,
491 const double mu);
493 bool GetIonMobility(const size_t ie, const size_t ib, const size_t ia,
494 double& mu) {
495 return GetEntry(ie, ib, ia, "IonMobility", m_iMob, mu);
496 }
497
499 bool SetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
500 const size_t ia, const double dl) {
501 return SetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
502 }
504 bool GetIonLongitudinalDiffusion(const size_t ie, const size_t ib,
505 const size_t ia, double& dl) {
506 return GetEntry(ie, ib, ia, "IonLongitudinalDiffusion", m_iDifL, dl);
507 }
509 bool SetIonTransverseDiffusion(const size_t ie, const size_t ib,
510 const size_t ia, const double dt) {
511 return SetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
512 }
514 bool GetIonTransverseDiffusion(const size_t ie, const size_t ib,
515 const size_t ia, double& dt) {
516 return GetEntry(ie, ib, ia, "IonTransverseDiffusion", m_iDifT, dt);
517 }
519 bool SetIonDissociation(const size_t ie, const size_t ib, const size_t ia,
520 const double diss) {
521 return SetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
522 }
524 bool GetIonDissociation(const size_t ie, const size_t ib, const size_t ia,
525 double& diss) {
526 return GetEntry(ie, ib, ia, "IonDissociation", m_iDis, diss);
527 }
528
530 bool SetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia,
531 const double mu);
533 bool GetNegativeIonMobility(const size_t ie, const size_t ib, const size_t ia,
534 double& mu) {
535 return GetEntry(ie, ib, ia, "NegativeIonMobility", m_nMob, mu);
536 }
537
539 virtual void ResetTables();
540
541 void ResetElectronVelocity() {
542 m_eVelE.clear();
543 m_eVelB.clear();
544 m_eVelX.clear();
545 m_eVelWv.clear();
546 m_eVelWr.clear();
547 }
548 void ResetElectronDiffusion() {
549 m_eDifL.clear();
550 m_eDifT.clear();
551 m_eDifM.clear();
552 }
553 void ResetElectronTownsend() { m_eAlp.clear(); }
554 void ResetElectronAttachment() { m_eAtt.clear(); }
555 void ResetElectronTOFRates() {
556 m_eRIon.clear();
557 m_eRAtt.clear();
558 }
559 void ResetElectronLorentzAngle() { m_eLor.clear(); }
560
561 void ResetHoleVelocity() {
562 m_hVelE.clear();
563 m_hVelB.clear();
564 m_hVelX.clear();
565 }
566 void ResetHoleDiffusion() {
567 m_hDifL.clear();
568 m_hDifT.clear();
569 m_hDifM.clear();
570 }
571 void ResetHoleTownsend() { m_hAlp.clear(); }
572 void ResetHoleAttachment() { m_hAtt.clear(); }
573
574 void ResetIonMobility() {
575 m_iMob.clear();
576 m_iVel.clear();
577 }
578 void ResetIonDiffusion() {
579 m_iDifL.clear();
580 m_iDifT.clear();
581 }
582 void ResetIonDissociation() { m_iDis.clear(); }
583 void ResetNegativeIonMobility() {
584 m_nMob.clear();
585 m_nVel.clear();
586 }
587
588 void VelocityFromMobility(
589 const std::vector<std::vector<std::vector<double> > >& mob,
590 std::vector<std::vector<std::vector<double> > >& vel);
591
594 void SetExtrapolationMethodVelocity(const std::string& extrLow,
595 const std::string& extrHigh);
596 void SetExtrapolationMethodDiffusion(const std::string& extrLow,
597 const std::string& extrHigh);
598 void SetExtrapolationMethodTownsend(const std::string& extrLow,
599 const std::string& extrHigh);
600 void SetExtrapolationMethodAttachment(const std::string& extrLow,
601 const std::string& extrHigh);
602 void SetExtrapolationMethodIonMobility(const std::string& extrLow,
603 const std::string& extrHigh);
604 void SetExtrapolationMethodIonDissociation(const std::string& extrLow,
605 const std::string& extrHigh);
606
608 void SetInterpolationMethodVelocity(const unsigned int intrp);
609 void SetInterpolationMethodDiffusion(const unsigned int intrp);
610 void SetInterpolationMethodTownsend(const unsigned int intrp);
611 void SetInterpolationMethodAttachment(const unsigned int intrp);
612 void SetInterpolationMethodIonMobility(const unsigned int intrp);
613 void SetInterpolationMethodIonDissociation(const unsigned int intrp);
614
615 // Scaling of fields and transport parameters.
616 virtual double ScaleElectricField(const double e) const { return e; }
617 virtual double UnScaleElectricField(const double e) const { return e; }
618 virtual double ScaleVelocity(const double v) const { return v; }
619 virtual double ScaleDiffusion(const double d) const { return d; }
620 virtual double ScaleDiffusionTensor(const double d) const { return d; }
621 virtual double ScaleTownsend(const double alpha) const { return alpha; }
622 virtual double ScaleAttachment(const double eta) const { return eta; }
623 virtual double ScaleLorentzAngle(const double lor) const { return lor; }
624 virtual double ScaleDissociation(const double diss) const { return diss; }
625
626 // Optical properties
628 virtual bool GetOpticalDataRange(double& emin, double& emax,
629 const unsigned int i = 0);
631 virtual bool GetDielectricFunction(const double e, double& eps1, double& eps2,
632 const unsigned int i = 0);
633 // Get the photoabsorption cross-section [cm2] at a given energy.
634 virtual bool GetPhotoAbsorptionCrossSection(const double e, double& sigma,
635 const unsigned int i = 0);
636 virtual double GetPhotonCollisionRate(const double e);
637 virtual bool PhotonCollision(const double e, int& type, int& level,
638 double& e1, double& ctheta,
639 std::vector<Secondary>& secondaries);
640
642 void EnableDebugging() { m_debug = true; }
643 void DisableDebugging() { m_debug = false; }
644
646 virtual double CreateGPUTransferObject(MediumGPU*& med_gpu);
647
648 protected:
649 std::string m_className = "Medium";
650
651 static int m_idCounter;
652
653 // Number of components
654 unsigned int m_nComponents = 1;
655 // Name
656 std::string m_name = "";
657 // Temperature [K]
658 double m_temperature = 293.15;
659 // Pressure [Torr]
660 double m_pressure = 760.;
661 // Static dielectric constant
662 double m_epsilon = 1.;
663 // (Effective) atomic number Z
664 double m_z = 1.;
665 // Atomic weight A
666 double m_a = 0.;
667 // Number density [cm-3]
668 double m_density = 0.;
669#endif
670
671 // Id number
672 int m_id;
673
674 // Transport flags
675 bool m_driftable = false;
676 bool m_microscopic = false;
677 bool m_ionisable = false;
678
679#ifndef __GPUCOMPILE__
680
681 // W value
682 double m_w = 0.;
683 // Fano factor
684 double m_fano = 0.;
685
686 // Update flag
687 bool m_isChanged = true;
688
689 // Switch on/off debugging messages
690 bool m_debug = false;
691
692 // Tables of transport parameters
693 bool m_tab2d = false;
694
695 // Field grids
696 std::vector<double> m_eFields;
697 std::vector<double> m_bFields;
698 std::vector<double> m_bAngles;
699
700 // Electrons
701 std::vector<std::vector<std::vector<double> > > m_eVelE;
702 std::vector<std::vector<std::vector<double> > > m_eVelX;
703 std::vector<std::vector<std::vector<double> > > m_eVelB;
704 std::vector<std::vector<std::vector<double> > > m_eDifL;
705 std::vector<std::vector<std::vector<double> > > m_eDifT;
706 std::vector<std::vector<std::vector<double> > > m_eAlp;
707 std::vector<std::vector<std::vector<double> > > m_eAtt;
708 std::vector<std::vector<std::vector<double> > > m_eLor;
709 std::vector<std::vector<std::vector<double> > > m_eVelWv;
710 std::vector<std::vector<std::vector<double> > > m_eVelWr;
711 std::vector<std::vector<std::vector<double> > > m_eRIon;
712 std::vector<std::vector<std::vector<double> > > m_eRAtt;
713
714 std::vector<std::vector<std::vector<std::vector<double> > > > m_eDifM;
715
716 // Holes
717 std::vector<std::vector<std::vector<double> > > m_hVelE;
718 std::vector<std::vector<std::vector<double> > > m_hVelX;
719 std::vector<std::vector<std::vector<double> > > m_hVelB;
720 std::vector<std::vector<std::vector<double> > > m_hDifL;
721 std::vector<std::vector<std::vector<double> > > m_hDifT;
722 std::vector<std::vector<std::vector<double> > > m_hAlp;
723 std::vector<std::vector<std::vector<double> > > m_hAtt;
724
725 std::vector<std::vector<std::vector<std::vector<double> > > > m_hDifM;
726
727 // Ions
728 std::vector<std::vector<std::vector<double> > > m_iMob;
729 std::vector<std::vector<std::vector<double> > > m_iVel;
730 std::vector<std::vector<std::vector<double> > > m_iDifL;
731 std::vector<std::vector<std::vector<double> > > m_iDifT;
732 std::vector<std::vector<std::vector<double> > > m_iDis;
733 // Negative ions
734 std::vector<std::vector<std::vector<double> > > m_nMob;
735 std::vector<std::vector<std::vector<double> > > m_nVel;
736
737 // Thresholds for Townsend, attachment and dissociation coefficients.
738 unsigned int m_eThrAlp = 0;
739 unsigned int m_eThrAtt = 0;
740 unsigned int m_hThrAlp = 0;
741 unsigned int m_hThrAtt = 0;
742 unsigned int m_iThrDis = 0;
743
744 // Extrapolation methods (TODO: enum).
745 std::pair<unsigned int, unsigned int> m_extrVel = {0, 1};
746 std::pair<unsigned int, unsigned int> m_extrDif = {0, 1};
747 std::pair<unsigned int, unsigned int> m_extrAlp = {0, 1};
748 std::pair<unsigned int, unsigned int> m_extrAtt = {0, 1};
749 std::pair<unsigned int, unsigned int> m_extrLor = {0, 1};
750 std::pair<unsigned int, unsigned int> m_extrMob = {0, 1};
751 std::pair<unsigned int, unsigned int> m_extrDis = {0, 1};
752
753 // Interpolation methods
754 unsigned int m_intpVel = 2;
755 unsigned int m_intpDif = 2;
756 unsigned int m_intpAlp = 2;
757 unsigned int m_intpAtt = 2;
758 unsigned int m_intpLor = 2;
759 unsigned int m_intpMob = 2;
760 unsigned int m_intpDis = 2;
761
762 bool Velocity(const double ex, const double ey, const double ez,
763 const double bx, const double by, const double bz,
764 const std::vector<std::vector<std::vector<double> > >& velE,
765 const std::vector<std::vector<std::vector<double> > >& velB,
766 const std::vector<std::vector<std::vector<double> > >& velX,
767 const double q, double& vx, double& vy, double& vz) const;
768 bool VelocityFluxBulk(
769 const double ex, const double ey, const double ez, const double bx,
770 const double by, const double bz,
771 const std::vector<std::vector<std::vector<double> > >& velWv,
772 const std::vector<std::vector<std::vector<double> > >& velWr,
773 double& wv, double& wr) const;
774 static void Langevin(const double ex, const double ey, const double ez,
775 double bx, double by, double bz, const double mu,
776 double& vx, double& vy, double& vz);
777 static void Langevin(const double ex, const double ey, const double ez,
778 double bx, double by, double bz, const double mu,
779 const double muH, double& vx, double& vy, double& vz);
780 bool Diffusion(const double ex, const double ey, const double ez,
781 const double bx, const double by, const double bz,
782 const std::vector<std::vector<std::vector<double> > >& difL,
783 const std::vector<std::vector<std::vector<double> > >& difT,
784 double& dl, double& dt) const;
785 bool Diffusion(
786 const double ex, const double ey, const double ez, const double bx,
787 const double by, const double bz,
788 const std::vector<std::vector<std::vector<std::vector<double> > > >& diff,
789 double cov[3][3]) const;
790 bool Alpha(const double ex, const double ey, const double ez, const double bx,
791 const double by, const double bz,
792 const std::vector<std::vector<std::vector<double> > >& tab,
793 unsigned int intp, const unsigned int thr,
794 const std::pair<unsigned int, unsigned int>& extr,
795 double& alpha) const;
796 double GetAngle(const double ex, const double ey, const double ez,
797 const double bx, const double by, const double bz,
798 const double e, const double b) const;
799 bool Interpolate(const double e, const double b, const double a,
800 const std::vector<std::vector<std::vector<double> > >& table,
801 double& y, const unsigned int intp,
802 const std::pair<unsigned int, unsigned int>& extr,
803 const bool logval = false) const;
804
805 double Interpolate1D(const double e, const std::vector<double>& table,
806 const std::vector<double>& fields,
807 const unsigned int intpMeth,
808 const std::pair<unsigned int, unsigned int>& extr,
809 const bool logval = false) const;
810
811 bool SetEntry(const size_t i, const size_t j, const size_t k,
812 const std::string& fcn,
813 std::vector<std::vector<std::vector<double> > >& tab,
814 const double val);
815 bool GetEntry(const size_t i, const size_t j, const size_t k,
816 const std::string& fcn,
817 const std::vector<std::vector<std::vector<double> > >& tab,
818 double& val) const;
819
820 void SetExtrapolationMethod(const std::string& low, const std::string& high,
821 std::pair<unsigned int, unsigned int>& extr,
822 const std::string& fcn);
823 bool GetExtrapolationIndex(std::string str, unsigned int& nb) const;
824 size_t SetThreshold(
825 const std::vector<std::vector<std::vector<double> > >& tab) const;
826
827 void Clone(std::vector<std::vector<std::vector<double> > >& tab,
828 const std::vector<double>& efields,
829 const std::vector<double>& bfields,
830 const std::vector<double>& angles, const unsigned int intp,
831 const std::pair<unsigned int, unsigned int>& extr,
832 const double init, const std::string& label);
833 void Clone(std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
834 const size_t n, const std::vector<double>& efields,
835 const std::vector<double>& bfields,
836 const std::vector<double>& angles, const unsigned int intp,
837 const std::pair<unsigned int, unsigned int>& extr,
838 const double init, const std::string& label);
839
840 void Init(const size_t nE, const size_t nB, const size_t nA,
841 std::vector<std::vector<std::vector<double> > >& tab,
842 const double val);
843 void Init(const size_t nE, const size_t nB, const size_t nA, const size_t nT,
844 std::vector<std::vector<std::vector<std::vector<double> > > >& tab,
845 const double val);
846
847#else
848
849#include "MediumMagboltz.hh"
850
851 friend class MediumGas;
852 friend class MediumMagboltz;
853
854 // enum to mimic polymorphism
855 enum class MediumType { Medium = 0, MediumGas, MediumMagboltz };
856
857 MediumType m_MediumType{MediumType::Medium};
858
859#endif
860};
#define __DEVICE__
Base class for gas media.
Definition MediumGas.hh:25
Interface to Magboltz (version 11).

◆ GARFIELD_CLASS_NAME() [3/4]

class Garfield::GARFIELD_CLASS_NAME ( Sensor )

Sensor

Default constructor.

Destructor.

Constructor from a single component.

Add a component.

Get the number of components attached to the sensor.

Retrieve the pointer to a given component.

Activate/deactivate a given component.

Activate/deactivate use of the magnetic field of a given component.

Does the sensor have a non-zero magnetic field?

Add an electrode.

Get the number of electrodes attached to the sensor.

Remove all electrodes.

Remove all components, electrodes and reset the sensor.

Get the drift field and potential at (x, y, z).

Get the drift field at (x, y, z).

Get the magnetic field at (x, y, z).

Get the weighting field at (x, y, z).

Get the weighting potential at (x, y, z).

Get the delayed weighting potential at (x, y, z).

Get the medium at (x, y, z).

Switch debugging messages on/off.

Set the user area to the default.

Set the user area explicitly.

Return the current user area.

Check if a point is inside the user area.

Check if a point is inside an active medium and inside the user area.

Return the voltage range.

Start a new event, when computing the average signal over multiple events.

Reset signals and induced charges of all electrodes.

Set the time window and binning for the signal calculation.

Parameters
tstartstart time [ns]
tstepbin width [ns]
nstepsnumber of bins

Retrieve the time window and binning.

Compute the component of the signal due to the delayed weighting field.

Set the points in time at which to evaluate the delayed weighting field.

Set the number of points to be used when averaging the delayed signal vector over a time bin (default: 0). The averaging is done with a $2\times navg + 1$ point Newton-Raphson integration.

Set/override the signal in a given time bin explicitly.

Set/override the signal.

Retrieve the total signal for a given electrode and time bin.

Retrieve the electron signal for a given electrode and time bin.

Retrieve the prompt signal for a given electrode and time bin.

Retrieve the delayed signal for a given electrode and time bin.

Retrieve the electron signal for a given electrode and time bin.

Retrieve the ion or hole signal for a given electrode and time bin.

Retrieve the delayed electron signal for a given electrode and time bin.

Retrieve the delayed ion/hole signal for a given electrode and time bin.

Calculated using the weighting potentials at the start and end points.

Set the function to be used for evaluating the transfer function.

Set the points to be used for interpolating the transfer function.

Set the transfer function using a Shaper object.

Evaluate the transfer function at a given time.

Print some information about the presently set transfer function.

Plot the presently set transfer function.

Cache integral and FFT of the transfer function instead of recomputing it at every call (default: on).

Convolute the induced current on a given electrode with the transfer function.

Convolute all induced currents with the transfer function.

Replace the signal on a given electrode by its integral.

Replace the signals on all electrodes curve by their integrals.

Return whether the signal has been integrated/convoluted.

Delay the signal and subtract an attenuated copy (modelling a constant fraction discriminator).

\[  v_{out} = v_{in}\left(t - t_d\right) - f v_{in}.
\]

Set the function to be used for evaluating the noise component.

Add noise to the induced signal.

Add white noise to the induced signal, given a desired output ENC.

Parameters
labelname of the electrode
encEquivalent Noise Charge, in electrons.
poissonflag to sample the number of noise pulses from a Poisson distribution. Otherwise the noise charge in each bin is sampled from a Gaussian distribution.
q0amplitude of the noise delta pulses (in electrons).

Add white noise to the induced signals on all electrodes.

Determine the threshold crossings of the current signal curve.

Parameters
thrthreshold value
labelelectrode for which to compute the threshold crossings
nnumber of threshold crossings

Get the number of threshold crossings (after having called ComputeThresholdCrossings).

Retrieve the time and type of a given threshold crossing (after having called ComputeThresholdCrossings.

Parameters
iindex
timethreshold crossing time [ns]
levelthreshold (should correspond to the value requested).
riseflag whether the crossing is on a rising or falling slope.

Calculate the signal from a drift line using the weighting potential method.

Calculate the signal from an avalanche using the weighting potential method.

Calculate the signal from a drift line using the weighting field method.

Calculate the signal from a drift line using the weighting field method, given the drift velocities at each drift line point.

Plot the induced signal.

Exporting induced signal to a csv file.

Add the induced charge from a charge carrier drift between two points.

Determine whether a line between two points crosses a wire, calls Component::CrossedWire.

Determine whether a point is in the trap radius of a wire.

Determine whether a line between two points crosses a plane, calls Component::CrossedPlane.

Integrate the electric field flux through a line from (x0,y0,z0) to (x1,y1,z1) along a direction (xp,yp,zp).

Create and initialise GPU Transfer class

Mutex.

Components

Electrodes

Definition at line 1 of file Sensor.hh.

44 {
45 public:
47 GARFIELD_CLASS_NAME(Sensor)() = default;
49 ~GARFIELD_CLASS_NAME(Sensor)() = default;
50
51#ifndef __GPUCOMPILE__
53 Sensor(Component* comp);
55 void AddComponent(Component* comp);
57 size_t GetNumberOfComponents() const { return m_components.size(); }
59 Component* GetComponent(const unsigned int i);
61 void EnableComponent(const unsigned int i, const bool on);
63 void EnableMagneticField(const unsigned int i, const bool on);
65 bool HasMagneticField() const;
66
68 void AddElectrode(Component* comp, const std::string& label);
70 size_t GetNumberOfElectrodes() const { return m_electrodes.size(); }
72 void ClearElectrodes();
74 void Clear();
75
77 void ElectricField(const double x, const double y, const double z, double& ex,
78 double& ey, double& ez, double& v, Medium*& medium,
79 int& status);
80#endif
83 void ElectricField(const double x, const double y, const double z, double& ex,
84 double& ey, double& ez,
85 GARFIELD_CLASS_NAME(Medium) * &medium,
86 int& status) __GPUCONST__;
87
88#ifndef __GPUCOMPILE__
90 void MagneticField(const double x, const double y, const double z, double& bx,
91 double& by, double& bz, int& status);
92
94 void WeightingField(const double x, const double y, const double z,
95 double& wx, double& wy, double& wz,
96 const std::string& label);
98 double WeightingPotential(const double x, const double y, const double z,
99 const std::string& label);
100
102 double DelayedWeightingPotential(const double x, const double y,
103 const double z, const double t,
104 const std::string& label);
105
107 Medium* GetMedium(const double x, const double y, const double z);
108
110 void EnableDebugging(const bool on = true) { m_debug = on; }
111
113 bool SetArea(const bool verbose = false);
115 bool SetArea(const double xmin, const double ymin, const double zmin,
116 const double xmax, const double ymax, const double zmax);
118 bool GetArea(double& xmin, double& ymin, double& zmin, double& xmax,
119 double& ymax, double& zmax);
120#endif
123 bool IsInArea(const double x, const double y, const double z) __GPUCONST__;
124
125#ifndef __GPUCOMPILE__
127 bool IsInside(const double x, const double y, const double z);
128
130 bool GetVoltageRange(double& vmin, double& vmax);
131
133 void NewSignal() { ++m_nEvents; }
135 void ClearSignal();
136#else
138 __device__ void AddSignal(const double q, const double t0, const double t1,
139 const double x0, const double y0, const double z0,
140 const double x1, const double y1, const double z1,
141 const bool integrateWeightingField,
142 const bool useWeightingPotential,
143 const int particle_idx);
144#endif
145
146#ifndef __GPUCOMPILE__
152 void SetTimeWindow(const double tstart, const double tstep,
153 const unsigned int nsteps);
155 void GetTimeWindow(double& tstart, double& tstep,
156 unsigned int& nsteps) const {
157 tstart = m_tStart;
158 tstep = m_tStep;
159 nsteps = m_nTimeBins;
160 }
161
163 void EnableDelayedSignal(const bool on = true) { m_delayedSignal = on; }
165 void SetDelayedSignalTimes(const std::vector<double>& ts);
170 void SetDelayedSignalAveragingOrder(const unsigned int navg) {
171 m_nAvgDelayedSignal = navg;
172 }
173
175 void SetSignal(const std::string& label, const unsigned int bin,
176 const double signal);
178 void SetSignal(const std::string& label, const std::vector<double>& ts,
179 const std::vector<double>& is);
181 double GetSignal(const std::string& label, const unsigned int bin);
183 double GetSignal(const std::string& label, const unsigned int bin,
184 const int comp);
186 double GetPromptSignal(const std::string& label, const unsigned int bin);
188 double GetDelayedSignal(const std::string& label, const unsigned int bin);
190 double GetElectronSignal(const std::string& label, const unsigned int bin);
192 double GetIonSignal(const std::string& label, const unsigned int bin);
194 double GetDelayedElectronSignal(const std::string& label,
195 const unsigned int bin);
197 double GetDelayedIonSignal(const std::string& label, const unsigned int bin);
199 double GetInducedCharge(const std::string& label);
200
202 void SetTransferFunction(std::function<double(double)>);
204 void SetTransferFunction(const std::vector<double>& times,
205 const std::vector<double>& values);
207 void SetTransferFunction(Shaper& shaper);
209 double GetTransferFunction(const double t);
211 void PrintTransferFunction();
213 void PlotTransferFunction();
216 void EnableTransferFunctionCache(const bool on = true) {
217 m_cacheTransferFunction = on;
218 }
221 bool ConvoluteSignal(const std::string& label, const bool fft = false);
223 bool ConvoluteSignals(const bool fft = false);
225 bool IntegrateSignal(const std::string& label);
227 bool IntegrateSignals();
229 bool IsIntegrated(const std::string& label) const;
230
236 bool DelayAndSubtractFraction(const double td, const double f);
238 void SetNoiseFunction(double (*f)(double t));
240 void AddNoise(const bool total = true, const bool electron = false,
241 const bool ion = false);
250 void AddWhiteNoise(const std::string& label, const double enc,
251 const bool poisson = true, const double q0 = 1.);
253 void AddWhiteNoise(const double enc, const bool poisson = true,
254 const double q0 = 1.);
255
261 bool ComputeThresholdCrossings(const double thr, const std::string& label,
262 int& n);
265 size_t GetNumberOfThresholdCrossings() const {
266 return m_thresholdCrossings.size();
267 }
275 bool GetThresholdCrossing(const unsigned int i, double& time, double& level,
276 bool& rise) const;
277
280 void AddSignalWeightingPotential(
281 const double q, const std::vector<double>& ts,
282 const std::vector<std::array<double, 3> >& xs);
285 void AddSignalWeightingPotential(
286 const double q, const std::vector<double>& ts,
287 const std::vector<std::array<double, 3> >& xs,
288 const std::vector<double>& qs);
291 void AddSignalWeightingField(const double q, const std::vector<double>& ts,
292 const std::vector<std::array<double, 3> >& xs,
293 const bool integrateWeightingField);
297 void AddSignalWeightingField(const double q, const std::vector<double>& ts,
298 const std::vector<std::array<double, 3> >& xs,
299 const std::vector<std::array<double, 3> >& vs,
300 const std::vector<double>& ns, const int navg);
301
303 void PlotSignal(const std::string& label, TPad* pad,
304 const std::string optTotal = "t",
305 const std::string optPrompt = "",
306 const std::string optDelayed = "");
308 void ExportSignal(const std::string& label, const std::string& filename,
309 const bool chargeCariers = false) const;
310
312 void AddInducedCharge(const double q, const double x0, const double y0,
313 const double z0, const double x1, const double y1,
314 const double z1);
315
318 bool CrossedWire(const double x0, const double y0, const double z0,
319 const double x1, const double y1, const double z1,
320 double& xc, double& yc, double& zc, const bool centre,
321 double& rc);
323 bool InTrapRadius(const double q0, const double x0, const double y0,
324 const double z0, double& xw, double& yw, double& rw);
327 bool CrossedPlane(const double x0, const double y0, const double z0,
328 const double x1, const double y1, const double z1,
329 double& xc, double& yc, double& zc);
330
333 double IntegrateFluxLine(const double x0, const double y0, const double z0,
334 const double x1, const double y1, const double z1,
335 const double xp, const double yp, const double zp,
336 const unsigned int nI, const int isign = 0);
337 // TODO!
338 double GetTotalInducedCharge(const std::string& label);
339
340 double StepSizeHint();
341
343 double CreateGPUTransferObject(SensorGPU*& sensor_gpu);
344#if USEGPU
346 void TransferGPUElectrodeSignals(SensorGPU*& sensor_gpu);
347#endif
348
349 private:
350 std::string m_className = "Sensor";
352 std::mutex m_mutex;
353#endif
354
356#ifdef __GPUCOMPILE__
357 ComponentGPU** m_components = nullptr;
358 size_t m_numComponents;
359 friend class Sensor;
360#else
361 std::vector<std::tuple<Component*, bool, bool> > m_components;
362#endif
363
364#ifndef __GPUCOMPILE__
365 struct Electrode {
366 Component* comp;
367 std::string label;
368 std::vector<double> signal;
369 std::vector<double> delayedSignal;
370 std::vector<double> electronSignal;
371 std::vector<double> ionSignal;
372 std::vector<double> delayedElectronSignal;
373 std::vector<double> delayedIonSignal;
374 double charge;
375 bool integrated;
376 };
378 std::vector<Electrode> m_electrodes;
379#else
380 struct ElectrodeGPU {
381 ComponentGPU* comp;
382 int label;
383 double* signal;
384 };
385
386 ElectrodeGPU* m_electrodes = nullptr;
387 size_t m_numElectrodes;
388#endif
389
390 // Time window for signals
391 double m_tStart = 0.;
392 double m_tStep = 10.;
393 unsigned int m_nTimeBins = 200;
394 unsigned int m_nEvents = 0;
395#ifndef __GPUCOMPILE__
396 bool m_delayedSignal = false;
397 std::vector<double> m_delayedSignalTimes;
398 unsigned int m_nAvgDelayedSignal = 0;
399
400 // Transfer function
401 std::function<double(double)> m_fTransfer;
402 Shaper* m_shaper = nullptr;
403 std::vector<std::pair<double, double> > m_fTransferTab;
404 bool m_cacheTransferFunction = true;
405 // Integral of the transfer function squared.
406 double m_fTransferSq = -1.;
407 // FFT of the transfer function.
408 std::vector<double> m_fTransferFFT;
409
410 // Noise
411 double (*m_fNoise)(double t) = nullptr;
412
413 std::vector<std::pair<double, bool> > m_thresholdCrossings;
414 double m_thresholdLevel = 0.;
415#endif
416
417 // User bounding box
418 double m_xMinUser = 0., m_yMinUser = 0., m_zMinUser = 0.;
419 double m_xMaxUser = 0., m_yMaxUser = 0., m_zMaxUser = 0.;
420 bool m_hasUserArea = false;
421
422#ifndef __GPUCOMPILE__
423 // Switch on/off debugging messages
424 bool m_debug = false;
425
426 // Return the current sensor size
427 bool GetBoundingBox(double& xmin, double& ymin, double& zmin, double& xmax,
428 double& ymax, double& zmax);
429
430 void FillSignal(Electrode& electrode, const double q,
431 const std::vector<double>& ts, const std::vector<double>& is,
432 const int navg, const bool delayed = false);
433 void FillBin(Electrode& electrode, const unsigned int bin,
434 const double signal, const bool electron, const bool delayed) {
435 std::lock_guard<std::mutex> guard(m_mutex);
436 electrode.signal[bin] += signal;
437 if (delayed) electrode.delayedSignal[bin] += signal;
438 if (electron) {
439 electrode.electronSignal[bin] += signal;
440 if (delayed) electrode.delayedElectronSignal[bin] += signal;
441 } else {
442 electrode.ionSignal[bin] += signal;
443 if (delayed) electrode.delayedIonSignal[bin] += signal;
444 }
445 }
446#else
447 __device__ void FillBin(ElectrodeGPU& electrode, const unsigned int bin,
448 const double signal, const bool electron,
449 const bool delayed, const int particle_idx);
450#endif
451
452#ifndef __GPUCOMPILE__
453 void IntegrateSignal(Electrode& electrode);
454 void ConvoluteSignal(Electrode& electrode, const std::vector<double>& tab);
455 bool ConvoluteSignalFFT();
456 bool ConvoluteSignalFFT(const std::string& label);
457 void ConvoluteSignalFFT(Electrode& electrode, const std::vector<double>& tab,
458 const unsigned int nn);
459 // Evaluate the integral over the transfer function squared.
460 double TransferFunctionSq();
461 double InterpolateTransferFunctionTable(const double t) const;
462 void MakeTransferFunctionTable(std::vector<double>& tab);
463 void FFT(std::vector<double>& data, const bool inverse, const int nn);
464#endif
465};
#define __GPUCONST__
Definition Sensor.hh:37
Class for signal processing.
Definition Shaper.hh:10

◆ GARFIELD_CLASS_NAME() [4/4]

class Garfield::GARFIELD_CLASS_NAME ( TetrahedralTree )

Helper class for searches in field maps.

This class stores the mesh nodes and elements in an Octree data structure to optimize the element search operations

Author: Ali Sheharyar

Organization: Texas A&M University at Qatar

Destructor

Insert a mesh element with given bounding box and index to the tree.

Create and initialise GPU Transfer class

Definition at line 1 of file TetrahedralTree.hh.

42 {
43 public:
44#ifdef __GPUCOMPILE__
45 // Constructor
46 GARFIELD_CLASS_NAME(TetrahedralTree)() = default;
47
48 // Destructor
49 ~GARFIELD_CLASS_NAME(TetrahedralTree)() = default;
50#else
51 // Constructor
52 GARFIELD_CLASS_NAME(TetrahedralTree)(const Vec3& origin,
53 const Vec3& halfDimension);
54
56 ~GARFIELD_CLASS_NAME(TetrahedralTree)();
57#endif
58
59#ifndef __GPUCOMPILE__
60 // Insert a mesh node (a vertex/point) to the tree
61 void InsertMeshNode(Vec3 point, const int index);
62
64 void InsertMeshElement(const double bb[6], const int index);
65
67 double CreateGPUTransferObject(TetrahedralTreeGPU*& tree_gpu);
68#endif
69
70 private:
71 static std::vector<int> emptyBlock;
72
73 // Physical centre of this tree node.
74 Vec3 m_origin;
75 // Half the width/height/depth of this tree node.
76 Vec3 m_halfDimension;
77 // Storing min and max points for convenience
78 Vec3 m_min, m_max;
79
80 // The tree has up to eight children and can additionally store
81 // a list of mesh nodes and mesh elements.
82 // Pointers to child octants.
83 GARFIELD_CLASS_NAME(TetrahedralTree) * children[8];
84
85 // Children follow a predictable pattern to make accesses simple.
86 // Here, - means less than 'origin' in that dimension, + means greater than.
87 // child: 0 1 2 3 4 5 6 7
88 // x: - - - - + + + +
89 // y: - - + + - - + +
90 // z: - + - + - + - +
91
92#ifndef __GPUCOMPILE__
93 std::vector<std::pair<Vec3, int> > nodes;
94#endif
95
96#ifdef __GPUCOMPILE__
97 int* elements{nullptr};
98 int numelements{0};
99#else
100 std::vector<int> elements;
101#endif
102
103 static const size_t BlockCapacity = 10;
104
105#ifndef __GPUCOMPILE__
106 // Check if the given box overlaps with this tree node.
107 bool DoesBoxOverlap(const double bb[6]) const;
108#endif
109 // Check if this tree node is a leaf or intermediate node.
110 __DEVICE__ bool IsLeafNode() const;
111
112 public:
113// Get all tetrahedra linked to a block corresponding to the given point
114#ifdef __GPUCOMPILE__
115 __device__ void GetElementsInBlock(const Vec3& point,
116 const int*& tet_list_elems,
117 int& num_elems) const;
118#else
119 const std::vector<int>& GetElementsInBlock(const Vec3& point) const;
120#endif
121
122 private:
123 __DEVICE__ int GetOctantContainingPoint(const Vec3& point) const;
124 // Get a block containing the input point
125 __DEVICE__ const GARFIELD_CLASS_NAME(TetrahedralTree) *
126 GetBlockFromPoint(const Vec3& point) const;
127 // A helper function used by the function above.
128 // Called recursively on the child nodes.
129 __DEVICE__ const GARFIELD_CLASS_NAME(TetrahedralTree) *
130 GetBlockFromPointHelper(const Vec3& point) const;
131
132#ifdef __GPUCOMPILE__
133 friend class TetrahedralTree;
134#endif
135};
Vec3Impl< double > Vec3

◆ ltrim()

void Garfield::ltrim ( std::string & line)
inline

Definition at line 11 of file Utilities.hh.

11 {
12 line.erase(line.begin(), std::find_if(line.begin(), line.end(), [](int ch) {
13 return !std::isspace(ch);
14 }));
15}

◆ RndmDirection()

void Garfield::RndmDirection ( double & dx,
double & dy,
double & dz,
const double length = 1. )
inline

Draw a random (isotropic) direction vector.

Definition at line 138 of file Random.hh.

139 {
140 const double phi = TwoPi * RndmUniform();
141 const double ctheta = 2 * RndmUniform() - 1.;
142 const double stheta = sqrt(1. - ctheta * ctheta);
143 dx = length * cos(phi) * stheta;
144 dy = length * sin(phi) * stheta;
145 dz = length * ctheta;
146}
double RndmUniform()
Draw a random number uniformly distributed in the range [0, 1).
Definition Random.hh:25

◆ RndmGaussian() [1/2]

double Garfield::RndmGaussian ( )
inline

Draw a Gaussian random variate with mean zero and standard deviation one.

Definition at line 35 of file Random.hh.

35 {
36 static bool cached = false;
37 static double u = 0.;
38 if (cached) {
39 cached = false;
40 return u;
41 }
42 // Box-Muller algorithm
43 u = 2. * RndmUniform() - 1.;
44 double v = 2. * RndmUniform() - 1.;
45 double r2 = u * u + v * v;
46 while (r2 > 1.) {
47 u = 2. * RndmUniform() - 1.;
48 v = 2. * RndmUniform() - 1.;
49 r2 = u * u + v * v;
50 }
51 const double p = sqrt(-2. * log(r2) / r2);
52 u *= p;
53 cached = true;
54 return v * p;
55}

◆ RndmGaussian() [2/2]

double Garfield::RndmGaussian ( const double mu,
const double sigma )
inline

Draw a Gaussian random variate with mean mu and standard deviation sigma.

Definition at line 58 of file Random.hh.

58 {
59 return mu + sigma * RndmGaussian();
60}
double RndmGaussian()
Draw a Gaussian random variate with mean zero and standard deviation one.
Definition Random.hh:35

◆ RndmHeedWF()

double Garfield::RndmHeedWF ( const double w,
const double f )

Draw a random energy needed to create a single electron in a material asymptotic work function W and Fano factor F, according to Igor Smirnov's phenomenological model.

◆ RndmLandau()

double Garfield::RndmLandau ( )

Draw a random number from a Landau distribution.

◆ RndmLorentzian()

double Garfield::RndmLorentzian ( const double mu,
const double gamma )
inline

Draw a Lorentzian random variate with mean mu and half-width at half maximum gamma.

Definition at line 64 of file Random.hh.

64 {
65 return mu + gamma * tan(Pi * (RndmUniform() - 0.5));
66}

◆ RndmPoisson()

int Garfield::RndmPoisson ( const double mean)

Draw a random number from a Poisson distribution.

◆ RndmPolya()

double Garfield::RndmPolya ( const double theta)
inline

Draw a Polya distributed random number.

Definition at line 101 of file Random.hh.

101 {
102 // Algorithm from Review of Particle Physics
103 // C. Amsler et al, Phys. Lett. B 667 (2008)
104 if (theta <= 0.) return -log(RndmUniformPos());
105 const double c = 3 * (theta + 1.) - 0.75;
106 double u1, u2, v1, v2, v3;
107 double x;
108 while (1) {
109 u1 = RndmUniformPos();
110 v1 = u1 * (1. - u1);
111 if (v1 == 0.) continue;
112 v2 = (u1 - 0.5) * sqrt(c / v1);
113 x = theta + v2;
114 if (x <= 0.) continue;
115 u2 = RndmUniformPos();
116 v3 = 64 * u2 * u2 * pow(v1, 3);
117 if (v3 <= 1. - 2 * v2 * v2 / x ||
118 log(v3) <= 2 * (theta * log(x / theta) - v2)) {
119 return x / (theta + 1.);
120 }
121 }
122}
double RndmUniformPos()
Draw a random number uniformly distributed in the range (0, 1).
Definition Random.hh:28

◆ RndmUniform()

double Garfield::RndmUniform ( )
inline

Draw a random number uniformly distributed in the range [0, 1).

Definition at line 25 of file Random.hh.

25{ return Random::Draw(); }
static double Draw() noexcept
Definition Random.hh:18

◆ RndmUniformPos()

double Garfield::RndmUniformPos ( )
inline

Draw a random number uniformly distributed in the range (0, 1).

Definition at line 28 of file Random.hh.

28 {
29 double r = RndmUniform();
30 while (r <= 0.) r = RndmUniform();
31 return r;
32}

◆ RndmVavilov()

double Garfield::RndmVavilov ( const double rkappa,
const double beta2 )

Draw a random number from a Vavilov distribution.

◆ RndmVoigt()

double Garfield::RndmVoigt ( const double mu,
const double sigma,
const double gamma )
inline

Draw a random number according to a Voigt function with mean mu.

The Voigt function is a convolution of a Gaussian (standard deviation sigma) and a Lorentzian (half width gamma).

Definition at line 72 of file Random.hh.

73 {
74 if (sigma <= 0.) return RndmLorentzian(mu, gamma);
75 const double a = gamma / (Sqrt2 * sigma);
76 const double x = RndmLorentzian(0., a) + RndmGaussian(0., 1. / Sqrt2);
77 return mu + x * Sqrt2 * sigma;
78}
double RndmLorentzian(const double mu, const double gamma)
Draw a Lorentzian random variate with mean mu and half-width at half maximum gamma.
Definition Random.hh:64

◆ RndmYuleFurry()

unsigned int Garfield::RndmYuleFurry ( const double mean)
inline

Draw a random number from a geometric distribution.

Definition at line 81 of file Random.hh.

81 {
82 if (mean <= 0.) return 0;
83 return 1 + static_cast<unsigned int>(std::log(RndmUniformPos()) /
84 std::log1p(-1. / mean));
85 /*
86 const double u = RndmUniform();
87 double p = 1. / mean;
88 const double q = 1. - p;
89 double sum = p;
90 unsigned int k = 1;
91 while (sum < u) {
92 ++k;
93 p *= q;
94 sum += p;
95 }
96 return k;
97 */
98}

◆ rtrim()

void Garfield::rtrim ( std::string & line)
inline

Definition at line 17 of file Utilities.hh.

17 {
18 line.erase(std::find_if(line.rbegin(), line.rend(),
19 [](int ch) { return !std::isspace(ch); })
20 .base(),
21 line.end());
22}

◆ SetDefaultStyle()

void Garfield::SetDefaultStyle ( )
inline

Definition at line 10 of file Plotting.hh.

10{ plottingEngine.SetDefaultStyle(); }
PlottingEngine plottingEngine

◆ SetSansSerif()

void Garfield::SetSansSerif ( )
inline

Definition at line 14 of file Plotting.hh.

14{ plottingEngine.SetSansSerif(); }

◆ SetSerif()

void Garfield::SetSerif ( )
inline

Definition at line 12 of file Plotting.hh.

12{ plottingEngine.SetSerif(); }

◆ startsWith()

bool Garfield::startsWith ( const std::string & line,
const std::string & s )
inline

Definition at line 33 of file Utilities.hh.

33 {
34 return (line.rfind(s, 0) == 0);
35}

◆ tokenize()

std::vector< std::string > Garfield::tokenize ( const std::string & line)
inline

Definition at line 24 of file Utilities.hh.

24 {
25 std::vector<std::string> words;
26 std::istringstream ss(line);
27 for (std::string word; ss >> word;) {
28 words.push_back(word);
29 }
30 return words;
31}

Variable Documentation

◆ gComponentNeBem3d

ComponentNeBem3d* Garfield::gComponentNeBem3d
extern

◆ plottingEngine

PlottingEngine Garfield::plottingEngine
extern