Garfield++ User Guide
Garfield++ User Guide
Garfield++ User Guide
Version 2017.1
H. Schindler
February 2017
Contents
1. Introduction 5
1.1. Class Structure . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 5
2. Getting Started 7
2.1. Installation . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 7
2.2. Examples . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 8
2.2.1. Drift Tube . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 9
2.2.2. GEM . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 11
3. Media 13
3.1. Transport Parameters . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 13
3.1.1. Transport Tables . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 14
3.2. Electron Scattering Rates . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 15
3.3. Gases . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 16
3.3.1. Ion Transport . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 17
3.3.2. Magboltz . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 18
3.4. Semiconductors . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 20
3.4.1. Transport Parameters . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 20
4. Components 23
4.1. Defining the Geometry . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 23
4.1.1. Visualizing the Geometry . . . . . . . . . . . . . . . . . . . . . . . . . . . 25
4.2. Field Maps . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 25
4.2.1. Ansys . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 25
4.2.2. Synopsys TCAD . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 27
4.2.3. Elmer . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 28
4.2.4. CST . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 28
4.2.5. Regular grids . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 29
4.2.6. Visualizing the Mesh . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 29
4.3. Analytic Fields . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 30
4.3.1. Describing the Cell . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 30
4.3.2. Periodicities . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 31
4.3.3. Cell Types . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 31
4.3.4. Weighting Fields . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 32
4.4. Other Components . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 33
4.5. Visualizing the Field . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 34
4.6. Sensor . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 35
5. Tracks 37
5.1. Heed . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 38
5.1.1. Delta Electron Transport . . . . . . . . . . . . . . . . . . . . . . . . . . . 39
5.1.2. Photon Transport . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 39
3
Contents 4
5.2. SRIM . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 39
6. Charge Transport 42
6.1. Runge-Kutta-Fehlberg Integration . . . . . . . . . . . . . . . . . . . . . . . . . . 42
6.2. Monte Carlo Integration . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 42
6.3. Microscopic Tracking . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 44
6.4. Visualizing Drift Lines . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 47
7. Signals 49
7.1. Readout Electronics . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . . 50
B. Gases 54
Bibliography 57
1. Introduction
Garfield++ is an object-oriented toolkit for the detailed simulation of particle detectors that use a
gas mixture or a semiconductor material as sensitive medium.
For calculating electric fields, three techniques are currently being offered:
• solutions in the thin-wire limit for devices made of wires and planes;
• interfaces with finite element programs, which can compute approximate fields in nearly
arbitrary two- and three-dimensional configurations with dielectrics and conductors;
• an interface with the Synopsys Sentaurus device simulation program [19].
In the future, an interface to the neBEM field solver [11, 12] (which already exists for Garfield
[21]), should be made available.
For calculating the transport properties of electrons in gas mixtures, an interface to the “Magboltz”
program [2, 3] is available.
The ionization pattern produced along the track of relativistic charged particles can be simulated
using the program “Heed” [18]. For simulating the ionization produced by low-energy ions, an
interface for importing results calculated using the SRIM software package [23] is available.
The present document aims to provide an overview of the Garfield++ classes and their key func-
tionalities, but does not provide an exhaustive description of all classes and functions. A number of
examples and code snippets are included which may serve as a basis for the user’s own programs.
Further examples and information can be found on the webpage http://cern.ch/garfieldpp. If you
have questions, doubts, comments etc. about the code or this manual, please don’t hesitate to
contact the authors. Any kind of constructive feedback is highly welcome.
An overview of the different types of classes is given in Fig. 1.1. Two main categories can be distin-
guished: classes for describing the detector (material properties, geometry, fields), and transport
classes which deal with tracing particles through the device. The two class types are linked by the
class Sensor.
The individual classes are explained in detail in the following chapters.
Readers familiar with the structure of (Fortran) Garfield [21] will recognize a rough correspondence
between the above classes and the sections of Garfield. Medium classes, for instance, can be
regarded as the counterpart of the &GAS section; Component classes are similar in scope to the
&CELL section.
Garfield++ also includes a number of classes for visualization purposes, e. g. for plotting drift lines,
making a contour plot of the electrostatic potential or inspecting the layout of the detector. These
classes rely extensively on the graphics classes of the ROOT framework [4].
5
Chapter 1. Introduction 6
Detector Description
material properties
• gas → Magboltz Medium Geometry
• silicon
field calculation
• analytic Component
• field maps
• neBEM
Sensor
charge transport
primary ionization
• microscopic
• Heed Track Drift • MC
• SRIM
Transport • RKF
Figure 1.1. Overview of the main classes in Garfield++ and their interplay.
2. Getting Started
2.1. Installation
The source code is hosted on a Subversion1 (svn) repository managed by the CERN Central SVN
service. Web interfaces for browsing the code are available at:
• http://svnweb.cern.ch/world/wsvn,
• http://svnweb.cern.ch/trac/garfield.
The following instructions describe how to download and build Garfield++ from source.
• Make sure that ROOT is installed. For installation instructions see https://root.cern.
ch/building-root or https://root.cern.ch/downloading-root.
• Define an environment variable GARFIELD_HOME pointing to the directory where the Garfield++
classes are to be located. In the following, we assume that we want to install Garfield in a
directory /home/mydir/garfield. If you are using bash, type
export GARFIELD_HOME=/home/mydir/garfield
Include the above lines also in the .bashrc (or .cshrc) file in your home directory. If unsure
which shell you are using, type echo $SHELL.
• Download (“check out”) the code from the repository. This can be done via SSH access or
via HTTP access. For SSH access, give the command
where username is your CERN afs login. For HTTPS access, give the command
1
For more information about Subversion, have a look at http://svn.web.cern.ch/svn/docs.php and the docu-
ments listed there.
7
Chapter 2. Getting Started 8
Alternatively, you can download the tarballs from the web interface (see the above address)
and extract them in the $GARFIELD_HOME directory.
• There are two options for building the library: (1) using directly the makefile in $GARFIELD_HOME,
or (2) using CMake.
– Alternatively, if you want to switch on debugging and switch off optimisation, type
export HEED_DATABASE=$GARFIELD_HOME/Heed/heed++/database/
svn update
followed by make (in case of trouble: try make clean; make), can be used for downloading the
latest version of the code from the repository.
2.2. Examples
Section 2.2.1 discusses the calculation of transport parameters with Magboltz, the use of ana-
lytic field calculation techniques, “macroscopic” simulation of electron and ion drift lines, and the
calculation of induced signals.
Microscopic transport of electrons and the use of finite element field maps are dealt with in
Sec. 2.2.2.
Sample macros and further examples can be found on the webpage (cern.ch/garfieldpp).
9 Chapter 2. Getting Started
Gas Table
First, we prepare a table of transport parameters (drift velocity, diffusion coefficients, Townsend
coefficient, and attachment coefficient) as a function of the electric field E (and, in general, also
the magnetic field B as well as the angle between E and B). In this example, we use a gas mixture
of 93% argon and 7% carbon dioxide at a pressure of 3 atm and room temperature.
We also have to specify the number of electric fields to be included in the table and the electric
field range to be covered. Here we use 20 field points between 100 V/cm and 100 kV/cm with
logarithmic spacing.
Now we run Magboltz to generate the gas table for this grid. As input parameter we have to specify
the number of collisions (in multiples of 107 ) over which the electron is traced by Magboltz.
This calculation will take a while, don’t panic. After the calculation is finished, we save the gas
table to a file for later use.
gas->WriteGasFile("ar_93_co2_7.gas");
Electric Field
For calculating the electric field inside the tube, we use the class ComponentAnalyticField which
can handle (two-dimensional) arrangements of wires, planes and tubes.
The Component requires a description of the geometry, that is a list of volumes and associated
media.
// Voltages
const double vWire = 3270.;
const double vTube = 0.;
// Add the wire in the center.
cmp->AddWire(0., 0., 2 * rWire, vWire, "s");
// Add the tube.
cmp->AddTube(rTube, vTube, 0, "t");
cmp->AddReadout("s");
we tell the Component to prepare the solution for the weighting field of the wire (which we have
given the label “s” before).
Finally we assemble a Sensor object which acts as an interface to the transport classes discussed
below.
In this (not very realistic) example, we want to calculate only the electron signal. We set the time
interval within which the signal is recorded by the sensor to 2 ns, with a binning of 0.02 ns.
Avalanche
For simulating the electron avalanche we use the class AvalancheMC which uses the previously
computed tables of transport parameters to calculate drift lines and multiplication.
11 Chapter 2. Getting Started
Using the class ViewSignal, we plot the current induced on the wire by the avalanche simulated
in the previous step.
2.2.2. GEM
Field Map
• loading the mesh (ELIST.lis, NLIST.lis), the list of nodal solutions (PRNSOL.lis), and
the material properties (MPLIST.lis);
• specifying the length unit of the values given in the .LIS files;
Gas
In this example, we calculate electron avalanches using “microscopic” Monte Carlo simulation, based
directly on the electron-atom/molecule cross-sections in the Magboltz database.
gas->SetMaxElectronEnergy(200.);
const bool verbose = true;
gas->Initialise(verbose);
In order to track a particle through the detector we have to tell ComponentAnsys123 which field
map material corresponds to which Medium.
Avalanche
For compound media (e. g. gas mixtures), the identifiers and fractions of the constituents are
available via
Medium classes provide functions for calculating the macroscopic transport parameters of electrons,
holes, and ions as a function of the electric and magnetic field:
bool ElectronVelocity(const double ex, const double ey, const double ez,
const double bx, const double by, const double bz,
double& vx, double& vy, double& vz);
bool ElectronDiffusion(const double ex, const double ey, const double ez,
const double bx, const double by, const double bz,
double& dl, double& dt);
bool ElectronTownsend(const double ex, const double ey, const double ez,
const double bx, const double by, const double bz,
double& alpha);
bool ElectronAttachment(const double ex, const double ey, const double ez,
const double bx, const double by, const double bz,
double& eta);
13
Chapter 3. Media 14
The above functions return true if the respective parameter is available at the requested field.
Analogous functions are available for holes (albeit of course not meaningful for gases), and also for
ions (except for the Townsend and attachment coefficients). A function specific to ions is
bool IonDissociation(const double ex, const double ey, const double ez,
const double bx, const double by, const double bz,
double& diss);
The transport parameters can either be stored in a one-dimensional table (as a function of the
electric field only) or in a three-dimensional table (as a function of E, B, and the angle θ between
E and B). If only a one-dimensional table is present and the drift velocity at B 6= 0 is requested,
the Laue-Langevin equation is used.
In the present version of the code, all transport parameters share the same grid of electric fields,
magnetic fields, and angles. By default, the field and angular ranges are
• 20 electric field points between 100 V / cm and 100 kV / cm, with logarithmic spacing
• B = 0, θ = 0
For specifying the field grid, two functions are available:
emin, emax min. and max. value of the electric field range to be covered by the tables
ne number of electric field grid points
logE flag specifying whether the E-field grid points should be evenly spaced (false), or logarith-
mically spaced (true)
bmin, bmax, ne magnetic field range and number of values
amin, amax, na angular range and number of angles
efields, bfields, angles lists of E, B, and θ (in ascending order)
Electric fields have to be supplied in V / cm, magnetic fields in Tesla, and angles in rad.
The gas tables are interpolated using Newton polynomials. The order of the interpolation polyno-
mials can be set by means of
For calculating electron drift lines using “microscopic tracking” (see Sec. 6.3), the preparation of
an electron transport table is not necessary, since this method is based directly on the electron-
Chapter 3. Media 16
The following functions which are meant to be called from within the class AvalancheMicroscopic
are available in Medium:
•
double GetElectronCollisionRate(const double e, const int band = 0);
returns the total scattering rate of an electron with energy e (in eV) in the Medium. The
band index is relevant only for semiconducting media.
•
bool GetElectronCollision(const double e, int& type, int& level,
double& e1, double& dx, double& dy, double& dz,
int& nion, int& ndx, int& band);
nion number of “ionisation products” (i. e. electrons and ions) created in the collision
3.3. Gases
There are currently two classes implemented which can be used for the description of gaseous
media: MediumGas and its daughter class MediumMagboltz. While MediumGas deals only with the
interpolation of gas tables and the import of gas files, MediumMagboltz – owing to an interface to
the Magboltz program [3] – allows the calculation of transport parameters. In addition, the latter
class provides access to the electron-molecule scattering cross-sections used in Magboltz and is
thus suitable for microscopic tracking (chapter 6).
A valid gas mixture comprises at least one and at most six different species. A list of the presently
available gases and their identifiers can be found in the appendix. The fractions have to be strictly
positive and may add up to any non-zero value; internally they will be normalized to one.
The gas density is specified in terms of pressure (Torr) and temperature (K):
Note that the density is calculated using the ideal gas law.
In the following example the gas mixture is set to Ar/CH4 (80/20) at atmospheric pressure and
20◦ C.
The function
void PrintGas();
prints information about the present transport parameter tables and cross-section terms (if avail-
able).
Extensive compilations of ion mobilities and diffusion coefficients can be found in Refs. [6–8,
22].
3.3.2. Magboltz
Magboltz, written by Steve Biagi, is a program [3] for the calculation of electron transport properties
in gas mixtures using semi-classical Monte Carlo simulation. It includes a database of electron-
atom/molecule cross-sections for a large number of detection gases.
MediumMagboltz allows running Magboltz for a given electric field, magnetic field and field an-
gle:
can be used.
Electron transport parameter tables can be saved to file and read from file by means of
The format of the gas file used in Garfield++ is compatible with the one used in Garfield 9.
Scattering Rates
As a prerequisite for “microscopic tracking” a table of the electron scattering rates (based on
the electron-atom/molecule cross-sections included in the Magboltz database) for the current gas
mixture and density needs to be prepared. This can be done using the function
If the flag verbose is set to true, some information (such as gas properties, and collision rates at
selected energies) is printed during the initialisation.
If
void EnableCrossSectionOutput();
is called prior to Initialise, a table of the cross-sections (as retrieved from Magboltz) is written
to a file cs.txt in the current working directory.
By default, the scattering rates table extends from 0 to 40 eV. The max. energy to be included in
the scattering rates table can be set using
int GetNumberOfLevels();
bool GetLevel(const int i, int& ngas, int& type, std::string& descr, double& e);
int GetNumberOfElectronCollisions();
int GetNumberOfElectronCollisions(int& nElastic, int& nIonising, int& nAttachment,
int& nInelastic, int& nExcitation, int& nSuperelastic);
Chapter 3. Media 20
The first function returns total number of electron collisions (i. e. calls to GetElectronCollisions)
since the last reset. The second function additionally provides the number of collisions of each
cross-section category (elastic, ionising etc.). The third function returns the number of collisions
for a specific cross-section term. The counters can be reset using
void ResetCollisionCounters();
Excitation Transfer
Penning transfer can be taken into account in terms of a transfer efficiency ri , i. e. the probability
for an excited level i with an excitation energy i exceeding the ionisation potential ion of the
mixture to be “converted” to an ionisation. The simulation of Penning transfer can be switched
on/off using
3.4. Semiconductors
Like for all Medium classes the user has the possibility to specify the transport parameters in
tabulated form as function of electric field, magnetic field, and angle. If no such tables have been
entered, the transport parameters are calculated based on empirical parameterizations (as used,
for instance, in device simulation programs). Several mobility models are implemented. For the
mobility µ0 at low electric fields, the following options are available:
• Using
21 Chapter 3. Media
electrons holes
µL [10−6 cm2 V−1 ns−1 ] β µL [10−6 cm2 V−1 ns−1 ] β
Sentaurus [10] 1.417 −2.5 0.4705 −2.5
Minimos [17] 1.43 −2.0 0.46 −2.18
Reggiani [14] 1.32 −2.0 0.46 −2.2
the values of low-field electron and hole mobilities can be specified explicitly by the user.
• The following functions select the model to be used for the mobility due to phonon scattering:
void SetLatticeMobilityModelMinimos();
void SetLatticeMobilityModelSentaurus();
void SetLatticeMobilityModelReggiani();
In all cases, the dependence of the lattice mobility µL on the temperature T is described by
β
T
µL (T ) = µL (T0 ) , T0 = 300 K.
T0
The values of the parameters µL (T0 ) and β used in the different models are shown in
Table 3.3. By default, the “Sentaurus” model is activated.
• The parameterization to be used for modelling the impact of doping on the mobility is specified
using
void SetDopingMobilityModelMinimos();
void SetDopingMobilityModelMasetti();
The first function activates the model used in Minimos 6.1 (see Ref. [17]). Using the second
function the model described in Ref. [13] is activated (default setting).
For modelling the velocity as function of the electric field, the following options are available:
• The method for calculating the high-field saturation velocities can be set using
The first function sets user-defined saturation velocities (in cm/ns) for electrons and holes.
The other functions activate different parameterizations for the saturation velocity as function
Chapter 3. Media 22
T0 0.87
e
vsat = 0.0107 cm/ns,
T
0.52
h T0
vsat = 0.00837 cm/ns,
T
where T0 = 300 K. The expressions for the other two implemented models can be found in
Refs. [14, 15].
• The parameterization of the mobility as function of the electric field to be used can be
selected using
void SetHighFieldMobilityModelMinimos();
void SetHighFieldMobilityModelCanali();
void SetHighFieldMobilityModelReggiani();
void SetHighFieldMobilityModelConstant();
The last function requests a constant mobility (i. e. linear dependence of the velocity on the
electric field). The models activated by the other functions used the following expressions
2µe0 µh0
µe (E) = e 2 , µh (E) = , (Minimos)
1 + vµh0
r
2µ E
1+ 1 + v e0 sat
sat
µe,h
µe,h (E) = 0
1 , (Canali [5])
e,h β e,h βe,h
µ0 E
1+ vsate,h
µe0 µh0
µe (E) = e 3/2 2/3 , µh (E) = h 2 1/2 , (Reggiani [14])
µ0 E µ E
1+ e
vsat 1 + v0h
sat
void SetImpactIonisationModelGrant();
void SetImpactIonisationModelVanOverstraetenDeMan();
The calculation of electric fields is done by classes derived from the abstract base class ComponentBase.
The key functions are
ex, ey, ez, v electric field and potential at the given position
m pointer to the medium at the given position; if there is no medium at this location, a null pointer
is returned
status status flag indicating where the point is located (see Table 4.1)
As mentioned above, the purpose of Component classes is to provide, for a given location, the
electric (and magnetic) field and a pointer to the Medium (if available). For the latter task, it is
obviously necessary to specify the geometry of the device. In case of field maps, the geometry is
already defined in the field solver. It is, therefore, sufficient to associate the materials of the field
map with the corresponding Medium classes.
For other components (e. g. analytic or user-parameterized fields), the geometry has to be defined
separately.
Simple structures can be described by the native geometry (GeometrySimple), which has only a
very restricted repertoire of shapes (solids). At present, the available solids are
value meaning
0 inside a drift medium
>0 inside a wire
-1 . . . -4 on the side of a plane where no wires are
-5 inside the mesh, but not in a drift medium
-6 outside the mesh
23
Chapter 4. Components 24
• SolidBox,
• SolidTube, and
• SolidSphere.
As an example, we consider a gas-filled tube with a diameter of 1 cm and a length of 20 cm (along
the z-axis), centred at the origin:
Solids may overlap. When the geometry is scanned (triggered, for instance, by calling GetMedium),
the first medium found is returned. The sequence of the scan is reversed with respect to the
assembly fo the geometry. Hence, the last medium added to the geometry is considered the
innermost.
For more complex structures, the class GeometryRoot can be used which provides an interface to
the ROOT geometry (TGeo). Using GeometryRoot, the above example would look like this:
In either case (GeometrySimple and GeometryRoot), after assembly the geometry is passed to
the Component as a pointer:
Similarly, the function ViewCell::Plot3d() paints a three-dimensional view of the cell layout.
4.2.1. Ansys
The interpolation of FEM field maps created with the program Ansys [1] is dealt with by the
classes ComponentAnsys121 and ComponentAnsys123. The class names refer to the type of mesh
element used by Ansys:
• ComponentAnys121 reads two-dimensional field maps with 8-node curved quadrilaterals (known
as “plane121” in Anys).
• ComponentAnsys123 reads three-dimensional field maps with quadric curved tetrahedra (known
as “solid123” in Ansys).
The field map is imported with the function
elist name of the file containing the list of elements (default: "ELIST.lis")
Chapter 4. Components 26
nlist name of the file containing the list of nodes (default: "NLIST.lis")
mplist name of the file containing the list of materials (default: "MPLIST.lis")
prnsol name of the file containing the nodal solutions (default: "PRNSOL.lis")
In order to enable charge transport and ionization, the materials in the map need to be associated
with Medium classes.
The materials in the field map are characterized by the relative dielectric constant ε and the
conductivity σ. These parameters are accessible through the functions
prnsol name of the file containing the nodal solution for the weighting field configuration
The weighting field map has to use the same mesh as the previously read “actual” field map.
For periodic structures, e. g. GEMs, one usually models only the basic cell of the geometry
and applies appropriate symmetry conditions to cover the whole problem domain. The available
symmetry operations are:
• simple periodicities,
• mirror periodicities,
• rotation symmetries.
Mirror periodicity and simple periodicity as well as axial periodicity and rotation symmetry are,
obviously, mutually exclusive. In case of axial periodicity, the field map has to cover an integral
fraction of 2π.
In order to assess the quality of the mesh, one can retrieve the dimensions of each mesh element
using
In the following example we make histograms of the aspect ratio and element size.
Electric fields calculated using the device simulation program Synopsys Sentaurus [19] can be
imported with the classes ComponentTcad2d and ComponentTcad3d.
datafilename name of the file containing the nodal solution; the filename typically typically ends
with _des.dat
Both files have to be exported in DF-ISE format, files in the default TDR format cannot be read.
To convert a TDR file to _.dat and .grid files, the Sentaurus tool tdx can be used
The classes have been tested with meshes created with the application Mesh which can pro-
duce axis-aligned two- and three-dimensional meshes. The only three-dimensional mesh elements
ComponentTcad3d can deal with are tetrahedra. A mesh which consists only of simplex elements
(triangles in 2D, tetrahedra in 3D), can be generated by invoking Mesh with the option -t.
After importing the files, the regions of the device where charge transport is to be enabled need to
be associated with Medium classes.
void EnablePeriodicityX();
void EnableMirrorPeriodicityX();
4.2.3. Elmer
The class ComponentElmer (contributed by J. Renner) allows one to import field maps created
with the open source field solver Elmer and the mesh tool Gmsh. A detailed tutorial can be found
on the webpage.
4.2.4. CST
The class ComponentCST (contributed by K. Zenker) reads field maps extracted from CST Studio.
More details can be found at http://www.desy.de/ zenker/FLC/garfieldpp.html.
29 Chapter 4. Components
Electric field values on a regular two-dimensional or three-dimensional grid can be imported using
the class ComponentVoxel. As a first step, the grid needs to be defined using the function
The electric field values (and potential) for each grid cell are read in using
withPotential flag whether the file contains an additional column with the electrostatic potential
withRegion flag whether the file contains an additional column with an integer value corresponding
to the region index (each region can be associated with a different medium)
scaleX, scalE, scaleP scaling factors to be applied to the coordinates, electric field values and
potentials
The available formats are XY, IJ, XYZ, and IJK, the first two for two-dimensional maps, and the
last two for three-dimensional maps. In case of XY (XYZ), the first two (three) columns contain
the x, y (and z) coordinates of a given point in the grid, followed by the electric field values (and
potential if available) at this point. The class then looks up the grid cell corresponding to this point
and assigns the electric field and potential accordingly. In case of IJ (IJK) the indices of the grid
cell along x, y (and z) are specified directly.
For visualizing the mesh imported from a FEM field map, the class ViewFEMesh (written by J.
Renner) is available. Using
bool ViewFEMesh::Plot();
Chapter 4. Components 30
then allows draws a two-dimensional projection of the drift lines stored in the ViewDrift class
together with the mesh. The plot can be customized using
colorid index of the ROOT color with which all elements of material matid are to be drawn
fill flag indicating whether to draw a wireframe mesh (false) or filled elements
As an illustration consider the following example (suppose that fm is a pointer to a field map
component and driftView is a pointer to a ViewDrift class)
For two-dimensional geometries consisting of wires, planes and tubes, semi-analytic calculation
techniques – based essentially on the capacitance matrix method – are implemented.
Wires, tubes and planes can be added to the cell layout by means of the following functions:
// Add a wire
void AddWire(const double x, const double y, const double d,
const double v, const std::string& label,
const double length = 100.,
const double tension = 50., const double rho = 19.3);
// Add a tube
void AddTube(const double r, const double v,
const int nEdges, const std::string& label);
// Add a plane at constant x
void AddPlaneX(const double x, const double v, const std::string& label);
// Add a plane at constant y
void AddPlaneY(const double y, const double v, const std::string& label);
31 Chapter 4. Components
In all of these functions, the potential v (in V) and a label (used for signal calculation) have to be
supplied as parameters.
For wires, the center of the wire (x, y) and its diameter (d) need to be specified. Optional
parameters are the wire length, the tension (more precisely, the mass in g of the weight used to
stretch the wire during the assembly) and the density (in g / cm3 ) of the wire material. These
parameters have no influence on the electric field. The number of wires that can be added is not
limited.
Tube-specific parameters are the radius1 (r) and the number of edges, which determines the shape
of the tube:
• n = 0: cylindrical pipe
• 3 ≤ n ≤ 8: regular polygon
There can be only one tube in a cell. The tube is always centered at the origin (0, 0).
Planes are described by their coordinates. A cell can have at most two x and two y planes. Planes
and tubes cannot be used together in the same cell.
The geometry can be reset (thereby deleting all wires, planes and tubes) by
void Clear();
Before assembling and inverting the capacitance matrix, a check is performed whether the provided
geometry matches the requirements. If necessary, the planes and wires are reordered. Wires outside
the tube or the planes as well as overlapping wires are removed.
4.3.2. Periodicities
The class supports simple periodicity in x and y direction. The periodic lengths are set by means
of
std::string GetCellType();
The weighting field calculation for a readout group – i. e. all elements (wires, planes, etc.) with
the same label – is activated by the function
In addition to the weighting fields of the elements used for the calculation of the (actual) electric
field, the weighting field for a strip segment of a plane can also be calculated. Strips can be defined
using
direction orientation of the strip (’y’ or ’z’ in case of x-planes, ’x’ or ’z’ in case of y -planes
x, y coordinate of the plane on which the strip is located
smin, smax min. and max. coordinate of the strip
The strip weighting field is calculated using an analytic expression for the field between two infinite
parallel plates which are kept at ground potential except for the strip segment, which is raised to
1 V. The anode-cathode distance d to be used for the evaluation of this expression can be set by
the user (variable gap in AddStripOnPlaneX, AddStripOnPlaneY). If this variable is not specified
(or set to a negative value), the following default values are used:
• if two planes are present in the cell, d is assumed to be the distance between those planes;
• if only one plane is present, d is taken to be the distance to the nearest wire.
Similarly, pixels can be defined using
Pixel weighting fields are calculated using the expressions given in Ref. [16].
For simple calculations, the class ComponentConstant can be used. As the name implies, it provides
a uniform electric field. The electric field and potential can be specified using
void SetElectricField(const double ex, const double ey, const double ez);
void SetPotential(const double x, const double p, const double z,
const double v);
void SetWeightingField(const double wx, const double wy, const double wz,
const std::string label);
The class ComponentUser takes the electric field and potential from a user-defined function.
// Depletion voltage
const double vdep = 160.;
// Applied voltage
const double v = 200.;
// Detector thickness
Chapter 4. Components 34
ey = ez = 0.;
ex = (v - vdep) / d + 2 * x * vdep / (d * d);
The class ViewField provides a number of functions for plotting the potential and field of a
component.
x0, . . . , z1 coordinates of the start and end points of the line along which the potential is to be
evaluated
option quantity to be plotted: potential ("v"/"p"/"phi"), magnitude of the electric field ("e"/"field"),
or individual components of the field ("ex", "ey", "ez")
drawopt option string passed on to the function Draw() of the ROOT TF2
The first two functions create a contour and surface plot, respectively, in the selected viewing
plane. The latter function plots the potential/field along the line (x0 , y0 , z0 ) → (x1 , y1 , z1 ).
The component or sensor of which the field is to be plotted is set by means of
The viewing plane and the region to be drawn can be specified using
n, nx, ny number of points in x and y direction (default for one-dimensional plots: n = 1000;
default for two-dimensional plots: nx = ny = 200)
The number of contour levels can be set using
For the voltage range, the (estimated) minimum and maximum values of the potential in the
component/sensor are used as defaults. The range of the electric and weighting fields should
always be set by the user.
4.6. Sensor
The Sensor class can be viewed as a composite of components. In order to obtain a complete
description of a detector, it is sometimes useful to combine fields from different Component classes.
For instance, one might wish to use a field map for the electric field, calculate the weighting
field using analytic methods, and use a parameterized B field. Superpositions of several electric,
magnetic and weighting fields are also possible.
Components are added using
While AddComponent tells the Sensor that the respective Component should be included in the
calculation of the electric and magnetic field, AddElectrode requests the weighting field named
label to be used for computing the corresponding signal.
To reset the sensor, thereby removing all components and electrodes, use
void Clear();
The total electric and magnetic fields (sum over all components) at a given position are accessible
through the functions ElectricField and MagneticField. The syntax is the same as for the
corresponding functions of the Component classes. Unlike the fields, materials cannot overlap. The
function Sensor::GetMedium, therefore, returns the first valid drift medium found.
Chapter 4. Components 36
xmin, . . . , zmax corners of the bounding box within which transport is enabled.
Calling SetArea() (without arguments) sets the user area to the envelope of all components (if
available).
In addition, the Sensor class takes care of signal calculations (Chapter 7).
5. Tracks
The purpose of classes of the type Track is to simulate ionization patterns produced by fast charged
particles traversing the detector.
Only particles which are sufficiently long lived to leave a track in a detector are considered. A list
of the available particles is given in Table 5.1.
The kinematics of the charged particle can be defined by means of a number of equivalent meth-
ods:
void NewTrack(const double x0, const double y0, const double z0,
const double t0,
const double dx0, const double dy0, const double dz0);
t0 starting time
37
Chapter 5. Tracks 38
The starting point of the track has to be inside an ionizable medium. Depending on the type of
Track class, there can be further restrictions on the type of Medium. If the specified direction
vector has zero norm, an isotropic random vector will be generated.
After successful initialization, the “clusters” produced along the track can be retrieved by
xcls, ycls, zcls, tcls position (and time) of the ionizing collision
The function returns false if the list of clusters is exhausted or if there is no valid track.
The concept of a “cluster” deserves some explanation. In the present context it refers to the
energy loss in a single ionizing collision of the primary charged particle and the secondary electrons
produced in this process.
5.1. Heed
The program Heed [18] is an implementation of the photo-absorption ionization (PAI) model. It
was written by I. Smirnov. An interface to Heed is available through the class TrackHeed.
After calling GetCluster, one can retrieve details about the electrons in the present cluster us-
ing
Heed simulates the energy degradation of δ electrons and the production of secondary (“conduc-
tion”) electrons using a phenomenological algorithm described in Ref. [18].
The asymptotic W value (eV) and the Fano factor of a Medium can be specified by the user by
means of the functions
If these parameters are not set, Heed uses internal default values. The default value for the Fano
factor is F = 0.19.
The transport of δ electrons can be activated or deactivated using
void EnableDeltaElectronTransport();
void DisableDeltaElectronTransport();
If δ electron transport is disabled, the number of electrons returned by GetCluster is the number
of “primary” ionisation electrons, i. e. the photo-electrons and Auger electrons. Their kinetic
energies and locations are accessible through the function GetElectron.
If δ electron transport is enabled (default setting), the function GetElectron returns the locations
of the “conduction” electrons as calculated by the internal δ transport algorithm of Heed. Since this
method does not provide the energy and direction of the secondary electrons, the corresponding
parameters in GetElectron are not meaningful in this case.
void TransportPhoton(const double x0, const double y0, const double z0,
const double t0, const double e0,
const double dx0, const double dy0, const double dz0,
int& nel);
5.2. SRIM
SRIM1 is a program for simulating the energy loss of ions in matter. It produces tables of stopping
powers, range and straggling parameters that can subsequently be imported in Garfield using the
1
Stopping and Range of Ions in Matter, www.srim.org
Chapter 5. Tracks 40
returns true if the SRIM output file was read successfully. The SRIM file contains the following
data
• a list of kinetic energies at which losses and straggling have been computed;
• average energy lost per unit distance via electromagnetic processes, for each energy;
• average energy lost per unit distance via nuclear processes, for each energy;
• projected path length, for each energy;
• longitudinal straggling, for each energy;
• transverse straggling, for each energy.
These can be visualized using the functions
void PlotEnergyLoss();
void PlotRange();
void PlotStraggling();
and printed out using the function TrackSrim::Print(). In addition to these tables, the file
also contains the mass and charge of the projectile, and the density of the target medium. These
properties are also imported and stored by TrackSrim when reading in the file. Unlike in case of
Heed, the particle therefore does not need to be specified by the user. In addition, the following
parameters need to be supplied "by hand" by the user.
• the work function (in eV) of the medium (TrackSrim::SetWorkFunction),
• the Fano factor of the medium (TrackSrim::SetFanoFactor),
• the atomic number Z and mass number A of the medium (TrackSrim::SetAtomicMassNumbers).
TrackSrim tries to generate individual tracks which statistically reproduce the average quantities
calculated by SRIM. Starting with the energy specified by the user, it iteratively
• computes (by interpolating in the tables) the electromagnetic and nuclear energy loss per
unit length at the current particle energy,
• calculates a step with a length over which the particle will produce on average a certain
number of electrons,
• updates the trajectory based on the longitudinal and transverse scatter at the current particle
energy,
• calculates a randomised actual electromagnetic energy loss over the step and updates the
particle energy.
This is repeated until the particle has no energy left or leaves the geometry. The model for
randomising the energy loss over a step can be set using the function
Model Description
0 No fluctuations
1 Untruncated Landau distribution
2 Vavilov distribution (provided the kinematic parameters are within the range of applicability,
otherwise fluctuations are disabled)
3 Gaussian distribution
4 Combination of Landau, Vavilov and Gaussian models,
each applied in their alleged domain of applicability
void EnableTransverseStraggling();
void DisableTransverseStraggling();
void EnableLongitudinalStraggling();
void DisableLongitudinalStraggling();
If energy loss fluctuations are used, longitudinal straggling should be disabled. By default, transverse
straggling is switched on and longitudinal straggling is switched off.
SRIM is aimed at low energy nuclear particles which deposit large numbers of electrons in a medium.
The grouping of electrons to a cluster is therefore somewhat arbitrary. By default, TrackSrim will
adjust the step size such that there are on average 100 clusters on the track. If the user specifies
a target cluster size, using
the step size will be chosen such that a cluster comprises on average n electrons. Alternatively, if
the user specifies a maximum number of clusters, using
the step size will be chosen such that on average there are n / 2 clusters on the track.
6. Charge Transport
On a phenomenological level, the drift of charge carriers under the influence of an electric field E
and a magnetic field B is described by the first order equation of motion
ṙ = vd E (r) , B (r) , (6.1)
where vd is the drift velocity. For the solution of (6.1), two methods are available in Garfield++:
For accurate simulations of electron trajectories in small-scale structures (with characteristic dimen-
sions comparable to the electron mean free path), and also for detailed calculations of ionisation
and excitation processes, transporting electrons on a microscopic level – i. e. based on the second-
order equation of motion – is the method of choice. Microscopic tracking of electrons is dealt with
by the class AvalancheMicroscopic.
• a step of length ∆s = vd ∆t in the direction of the drift velocity vd at the local field is
calculated (with either the time step ∆t or the distance ∆s being specified by the user);
• a random diffusion step is√sampled from three uncorrelated Gaussian distributions with stan-
dard deviation σL =√ DL ∆s for the component parallel to the drift velocity and standard
deviation σT = DT ∆s for the two transverse components;
• the two steps are added vectorially and the location is updated.
In the first case the integration is done using fixed time steps (default: 20 ps), in the second case
using fixed distance steps (default: 10 µm). Calling the third function instructs the class to do the
42
43 Chapter 6. Charge Transport
integration with exponentially distributed time steps with a mean equal to the specified multiple of
the “collision time”
mvd
τ= .
qE
The third method is activated by default.
bool DriftElectron(const double x0, const double y0, const double z0,
const double t0);
bool DriftHole(const double x0, const double y0, const double z0,
const double t0);
bool DriftIon(const double x0, const double y0, const double z0,
const double t0);
bool AvalancheElectron(const double x0, const double y0, const double z0,
const double t0, const bool hole = false);
bool AvalancheHole(const double x0, const double y0, const double z0,
const double t0, const bool electron = false);
bool AvalancheElectronHole(const double x0, const double y0, const double z0,
const double t0);
The flags hole and electron specify whether multiplication of holes and electrons, respec-
tively, should be taken into account in the simulation. In case of gas-based detectors, only
AvalancheElectron with hole = false is meaningful.
The starting and endpoints of electrons in the avalanche can be retrieved using
int GetNumberOfElectronEndpoints();
void GetElectronEndpoint(const int i,
double& x0, double& y0, double& z0, double& t0,
double& x1, double& y1, double& z1, double& t1,
int& status) const;
status status code indicating why the tracking of the electron was stopped.
Chapter 6. Charge Transport 44
void EnableMagneticField();
instructs the class to consider not only the electric but also the magnetic field in the evaluation of
the transport parameters. By default, magnetic fields are not taken into account.
For debugging purposes, attachment and diffusion can be switched off using
void DisableAttachment();
void DisableDiffusion();
void UnsetTimeWindow();
Microscopic tracking is (at present) only possible for electrons. It is implemented in the class
AvalancheMicroscopic. A calculation is started by means of
void AvalancheElectron(const double x0, const double y0, const double z0,
const double t0, const double e0,
const double dx0 = 0., const double dy0 = 0., const double dz0 = 0.);
int GetNumberOfElectronEndpoints();
void GetElectronEndpoint(const int i,
double& x0, double& y0, double& z0, double& t0, double& e0,
double& x1, double& y1, double& z1, double& t1, double& e1,
int& status);
x0, y0, z0, t0, e0 initial position, time and energy of the electron
x1, y1, z1, t1, e1 final position, time and energy of the electron
status status code indicating why the tracking of the electron was stopped.
The function
bool DriftElectron(const double x0, const double y0, const double z0,
const double t0, const double e0,
const double dx0, const double dy0, const double dz0);
traces only the initial electron but not the secondaries produced along its drift path (the input
parameters are the same as for AvalancheElectron).
After each collision, the histogram is filled with the current electron energy.
The effect of magnetic fields can be included in the stepping algorithm using the function
void EnableMagneticField();
Chapter 6. Charge Transport 46
By default, magnetic fields are not taken into account in the calculation.
Using
the size of an electron avalanche can be limited. After the avalanche has reached the specified
max. size, no further secondaries are added to the stack of electrons to be transported.
Like in AvalancheMC a time window can be set/unset using
The function specified in SetUserHandleStep is called prior to each free-flight step. The param-
eters passed to this function are
x, y, z, t position and time,
e energy before the step
dx, dy, dz direction,
47 Chapter 6. Charge Transport
For plotting drift lines and tracks the class ViewDrift can be used. After attaching a ViewDrift
object to a transport class, e. g. using
ViewDrift stores the trajectories which are calculated by the transport class. The drawing of the
trajectories is triggered by the function
void ViewDrift::Plot();
Chapter 6. Charge Transport 48
In case of AvalancheMicroscopic, it is usually not advisable to plot every single collision. The
number of collisions to be skipped for plotting can be set using
Signals are calculated using the Shockley-Ramo theorem. The current i (t) induced by a particle
with charge q at a position r moving at a velocity v is given by
where Ew is the so-called weighting field for the electrode to be read out.
The basic steps for calculating the current induced by the drift of electrons and ions/holes are:
1. Prepare the weighting field for the electrode to be read out. This step depends on the field
calculation technique (i. e. the type of Component) which is used (see Chapter 4).
2. Tell the Sensor that you want to use this weighting field for the signal calculation.
where cmp is a pointer to the Component which calculates the weighting field, and label (in
our example "readout") is the name you have assigned to the weighting field in the previous
step.
3. Setup the binning for the signal calculation.
The first parameter in this function is the lower time limit (in ns), the second one is the bin
width (in ns), and the last one is the number of time bins.
4. Switch on signal calculation in the transport classes using
void AvalancheMicroscopic::EnableSignalCalculation();
void AvalancheMC::EnableSignalCalculation();
The Sensor then records and accumulates the signals of all avalanches and drift lines which
are simulated.
5. The calculated signal can be retrieved using
The functions GetElectronSignal and GetIonSignal return the signal induced by negative
and positive charges, respectively. GetSignal returns the sum of both electron and hole
signals.
49
Chapter 7. Signals 50
void Sensor::ClearSignal();
For plotting the signal, the class ViewSignal can be used. As an illustration of the above recipe
consider the following example.
// Electrode label
const std::string label = "readout";
// Setup the weighting field.
// In this example we use a FEM field map.
ComponentAnsys123* fm = new ComponentAnsys123();
...
fm->SetWeightingField("WPOT.lis", label);
In order to model the signal-processing by the front-end electronics, the “raw signal” – i. e. the
induced current – can be convoluted with a so-called “transfer function”. The transfer function to
be applied can be set using
in which case the transfer function will be calculated by interpolation of the values provided in the
table.
In order to convolute the presently stored signal with the transfer function (specified using the
above function), the function
bool Sensor::ConvoluteSignal();
can be called.
As an example, consider the following transfer function
t
f (t) = e e1−t/τ , τ = 25 ns
τ
double transfer(double t) {
The basic units are cm for distances, g for (macroscopic) masses, and ns for times. Particle
energies, momenta and masses are expressed in eV, eV/c and eV/c 2 , respectively. For example,
the electron mass is given in eV/c 2 , whereas the mass density of a material is given in g/cm3 . The
mass of an atom is specified in terms of the atomic mass number A.
There are a few exceptions from this system of units, though.
• The unit for the magnetic field B corresponding to the above system of units (10−5 Tesla)
is impractical. Instead, magnetic fields are expressed in Tesla.
• Pressures are specified in Torr.
• Electric charge is expressed in fC.
A summary of commonly used quantities and their units is given in Table A.1.
The values of the physical constants used in the code are defined in the file FundamentalConstants.hh.
52
53 Appendix A. Units and Constants
Table B.1 shows a list of the gases available in the current version of Magboltz. The star rating
represents an estimate of the reliability of the cross-section data for the respective gas. A rating of
“5*” indicates a detailed, well-validated description of the cross-sections, while “2*” indicates a low
quality, that is a coarse modelling of the cross-sections associated with large uncertainties.
54
55 Appendix B. Gases
cC3 H6 cyclopropane 4*
CH3 OH methanol 3*
C2 H5 OH ethanol 3*
C3 H7 OH isopropanol 3*
C 3 H8 O 2 methylal 2*
C4 H10 O2 DME 4*
CF4 tetrafluoromethane 5*
CHF3 fluoroform 3*
C2 F 6 hexafluoroethane 4*
C 2 H2 F 4 tetrafluoroethane 2*
C3 F 8 octafluoropropane 3*
SF6 sulfur hexafluoride 3*
BF3 boron trifluoride 4*
CF3 Br bromotrifluoromethane 3*
NH3 ammonia 4*
N(CH3 )3 TMA 3*
SiH4 silane 4*
GeH4 germane 3*
[3] , Monte Carlo simulation of electron drift and diffusion in counting gases under the
influence of electric and magnetic fields, Nucl. Instr. Meth. A 421 (1999), 234–240.
[4] R. Brun, F. Rademakers, et al., ROOT: An Object-Oriented Data Analysis Framework, http:
//root.cern.ch.
[5] C. Canali, G. Majni, R. Minder, and G. Ottaviani, Electron and Hole Drift VelocityMeasure-
ments in Silicon and Their Empirical Relation to Electric Field and Temperature, IEEE Trans.
Electron Devices 22 (1975), 1045–1047.
[6] H. W. Ellis et al., Transport properties of gaseous ions over a wide energy range, At. Data
Nucl. Data Tables 17 (1976), 177–210.
[7] , Transport properties of gaseous ions over a wide energy range. Part II, At. Data
Nucl. Data Tables 22 (1978), 179–217.
[8] , Transport properties of gaseous ions over a wide energy range. Part III, At. Data
Nucl. Data Tables 31 (1984), 113–151.
[9] W. N. Grant, Electron and Hole Ionization Rates in Epitaxial Silicon at High Fields, Solid
State Electronics 16 (1973), 1189–1203.
[10] C. Lombardi, S. Manzini, A. Saporito, and M. Vanzi, A Physically Based Mobility Model for
Numerical Simulation of Nonplanar Devices, IEEE Trans. CAD 7 (1988), 1164–1171.
[13] G. Masetti, M. Severi, and S. Solmi, Modeling of Carrier Mobility Against Carrier Concen-
tration in Arsenic-, Phosphorus-, and Boron-Doped Silicon, IEEE Trans. Electron Devices 30
(1983), 764–769.
[14] M. A. Omar and L. Reggiani, Drift and diffusion of charge carriers in silicon and their empirical
relation to the electric field, Solid State Electronics 30 (1987), 693–697.
[16] W. Riegler and G. Aglieri Rinella, Point charge potential and weighting field of a pixel or pad
in a plane condenser, Nucl. Instr. Meth. A 767 (2014), 267 – 270.
57
Bibliography 58
[17] S. Selberherr, W. Hänsch, M. Seavey, and J. Slotboom, The Evolution of the MINIMOS
Mobility Model, Solid State Electronics 33 (1990), 1425–1436.
[18] I. B. Smirnov, Modeling of ionization produced by fast charged particles in gases, Nucl. Instr.
Meth. A 554 (2005), 474–493.
[20] R. van Overstraeten and H. de Man, Measurement of the Ionization Rates in Diffused Silicon
p-n Junctions, Solid State Electronics 13 (1970), 583–608.
[22] L. Viehland and E. A. Mason, Transport properties of gaseous ions over a wide energy range,
IV, At. Data Nucl. Data Tables 60 (1995), 37–95.
[23] J. F. Ziegler, J. Biersack, and U. Littmark, The Stopping and Range of Ions in Matter,
Pergamon Press, 1985.