SolidStateDetectors.SSD_examplesConstant
SSD_examples::Dict{Symbol,String}

Dictionary with the paths to the example detector configuration files provided by the package.

Find the possible keys of the dictionary with keys(SSD_examples).

The example detector configuration files can be loaded via

path_to_config_file = SSD_examples[:InvertedCoax]
sim = Simulation(path_to_config_file)
SolidStateDetectors.ADLChargeDriftModelType
ADLChargeDriftModel{T <: SSDFloat, M <: AbstractDriftMaterial, N, TM <: AbstractTemperatureModel{T}} <: AbstractChargeDriftModel{T}

Charge drift model for electrons and holes based on the AGATA Detector Library. Find a detailed description of the calculations in ADL Charge Drift Model.

Fields

  • electrons::CarrierParameters{T}: Parameters to describe the electron drift along the <100> and <111> axes.
  • holes::CarrierParameters{T}: Parameters to describe the hole drift along the <100> and <111> axes.
  • crystal_orientation::SMatrix{3,3,T,9}: Rotation matrix that transforms the global coordinate system to the crystal coordinate system given by the <100>, <010> and <001> axes of the crystal.
  • γ::SVector{N,SMatrix{3,3,T,9}}: Reciprocal mass tensors to the N valleys of the conduction band.
  • parameters::ADLParameters{T}: Parameters needed for the calculation of the electron drift velocity.
  • temperaturemodel::TM: Models to scale the resulting drift velocities with respect to temperature

See also CarrierParameters.

SolidStateDetectors.AbstractImpurityDensityType
abstract type AbstractImpurityDensity{T <: SSDFloat} end

Struct defining an impurity density inside a Semiconductor.

For each impurity density, there should be a method get_impurity_density which returns the impurity density in SI units (1/m³) at a given point.

Examples

Note

The sign of the impurity density is important. It is taken into account in the conversion to a charge density and, thus, defines where the semiconductor is n-type or p-type.

SolidStateDetectors.CarrierParametersType
struct CarrierParameters{T <: SSDFloat}

Parameters needed to describe the electron or hole drift along the <100> and <111> axes.

Fields

  • axis100::VelocityParameters{T}: Parameters to describe the charge drift along the <100> axis.
  • axis111::VelocityParameters{T}: Parameters to describe the charge drift along the <111> axis.

See also VelocityParameters.

SolidStateDetectors.ConstantChargeDensityType
struct ConstantChargeDensity{T <: SSDFloat} <: AbstractChargeDensity{T}

Charge density model that assumes a constant charge density everywhere.

Fields

  • ρ::T: the constant value of the charge density.

Definition in Configuration File

A ConstantChargeDensity is defined in the configuration file through the field charge_density (of a passive or surrounding) with name: constant and a value for ρ.

An example definition of a constant charge density looks like this:

charge_density:
  name: constant
  value: 1.0e-10 # C/m³
SolidStateDetectors.ConstantImpurityDensityType
struct ConstantImpurityDensity{T <: SSDFloat} <: AbstractImpurityDensity{T}

Impurity density model that assumes a constant impurity density everywhere.

Fields

  • ρ::T: the constant value of the impurity density.

Definition in Configuration File

A ConstantImpurityDensity is defined in the configuration file through the field impurity_density (of a semiconductor) with name: constant and a value for ρ.

An example definition of a constant impurity density looks like this:

impurity_density:
  name: constant
  value: 1.0e10 # 1/m³
SolidStateDetectors.ContactType
struct Contact{T, G, MT} <: AbstractContact{T}

Contact of a SolidStateDetector.

For the simulation of the ElectricPotential, all contacts are fixed to a constant potential value.

Parametric types:

  • T: Precision type.
  • G: Type of geometry.
  • MT: Type of material.

Fields

  • potential::T: Potential (in V) to which the contact will be fixed during the calculation of the ElectricPotential.
  • material::MT: Material of the contact.
  • id::Int: Unique id that will unambiguously identify the contact.
  • name::String: Custom name for the contact, relevant for plotting.
  • geometry::G: Geometry of the contact, see Constructive Solid Geometry (CSG).

Definition in Configuration File

A Contact is defined in the configuration file through an entry in the contacts array of a detector. It needs id, potential and geometry and can optionally be given a name and material.

An example definition of contacts looks like this:

contacts:
  - name: "n+ contact"
    id: 1
    potential: 5000V
    material: HPGe # optional
    geometry: # ....
  - name: "p+ contact"
    id: 2
    potential: 0
    material: HPGe #optional
    geometry: # ....
SolidStateDetectors.CylindricalChargeDensityType
struct CylindricalChargeDensity{T <: SSDFloat} <: AbstractChargeDensity{T}

Charge density model which assumes a linear gradient in charge density in each spatial dimension of a cylindrical coordinate system.

Fields

  • offsets::NTuple{3,T}: charge density values at the origin of each dimension.
  • gradients::NTuple{3,T}: linear slopes in r and z direction.

Definition in Configuration File

A CylindricalChargeDensity is defined in the configuration file through the field charge_density (of a passive or surrounding) with name: cylindrical and optional fields r and z that can each contain init for initial values at 0 and gradient for gradients in that dimension.

An example definition of a cylindrical charge density looks like this:

charge_density:
  name: cylindrical
  r:  # impurity profile with linear gradient in r
    init: 1.0e-10     # C/m³
    gradient: 1.0e-11 # C/m⁴
SolidStateDetectors.CylindricalImpurityDensityType
struct CylindricalImpurityDensity{T <: SSDFloat} <: AbstractImpurityDensity{T}

Impurity density model which assumes a linear gradient in impurity density in each spatial dimension of a cylindrical coordinate system.

Fields

  • offsets::NTuple{3,T}: impurity density values at the origin of each dimension.
  • gradients::NTuple{3,T}: linear slopes in r and z direction.

Definition in Configuration File

A CylindricalImpurityDensity is defined in the configuration file through the field impurity_density (of a passive or surrounding) with name: cylindrical and optional fields r and z that can each contain init for initial values at 0 and gradient for gradients in that dimension.

An example definition of a cylindrical impurity density looks like this:

impurity_density:
  name: cylindrical
  r:  # impurity profile with linear gradient in r
    init: 1.0e10     # 1/m³
    gradient: 1.0e11 # 1/m⁴
SolidStateDetectors.DeadVolumeType
struct DeadVolume{T, G} <: AbstractVirtualVolume{T}

Volume inside which the charge drift is set to zero.

Parametric types

  • T: Precision type.
  • G: Type of geometry.

Fields

Definition in Configuration File

A DeadVolume is defined through an entry in the virtual_drift_volumes array of a detector with model: dead. It needs a geometry and can optionally be given a name.

An example definition of dead volumes looks like this:

virtual_drift_volume:
  - name: Volume 1
    model: dead
    geometry: # ...
  - name: Volume 2
    model: dead
    geometry: # ...
SolidStateDetectors.DielectricDistributionType
struct DielectricDistribution{T, N, S, AT} <: AbstractArray{T, N}

Dielectric distribution, or distribution of the relative permittivity, $\epsilon_r$, needed to calculate the ElectricPotential. The dielectric distribution is unitless.

Parametric types

  • T: Element type of data.
  • N: Dimension of the grid and data array.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • data::Array{T, N}: Array containing the values of the dielectric distribution at the discrete points of the grid.
  • grid::Grid{T, N, S, AT}: Grid defining the discrete points at which the electric potential is determined.
Note

The data array contains the values of the dielectric distribution at the discrete points between the points defined by the axes ticks of the extended grid of grid.

Thus, size(data) == size(grid) .+ 1 !

SolidStateDetectors.DiscreteAxisType
struct DiscreteAxis{T, BL, BR, I} <: AbstractAxis{T, BL, BR, I}

Axis with discrete ticks which is used to define a dimension of a Grid.

Parametric types

  • T: Type of ticks
  • BL: Boundary condition at the left endpoint.
  • BR: Boundary condition at the right endpoint.
  • I: IntervalSets.Interval (closed or open boundaries)

The boundary conditions of a DiscreteAxis can be BL, BR ∈ {:periodic, :reflecting, :infinite, :r0, :fixed}.

Fields

  • interval::I: Interval that defines the range of the axis.
  • ticks::Vector{T}: Array of values that correspond to the discrete ticks of the axis.

See also Grid.

SolidStateDetectors.DiscreteAxisMethod
DiscreteAxis(left_endpoint::T, right_endpoint::T, BL::Symbol, BR::Symbol, L::Symbol, R::Symbol, ticks::AbstractVector{T}) where {T}

Constructor of a DiscreteAxis.

Arguments

  • left_endpoint::T: Left endpoint of the interval of the DiscreteAxis.
  • right_endpoint::T: Right endpoint of the interval of the DiscreteAxis.
  • BL::Symbol: Boundary condition at the left endpoint.
  • BR::Symbol: Boundary condition at the right endpoint.
  • L::Symbol: Boundary type of the left endpoint.
  • R::Symbol: Boundary type of the right endpoint.
  • ticks::AbstractVector{T}: Array of values that correspond to the discrete ticks of the axis.

The boundary conditions of a DiscreteAxis can be BL, BR ∈ {:periodic, :reflecting, :infinite, :r0, :fixed}.

The boundary types of a DiscreteAxis can be L, R ∈ {:closed, :open}.

Examples

DiscreteAxis(-2.0, 2.0, :infinite, :infinite, :closed, :closed, collect(-2:0.1:2))
SolidStateDetectors.EffectiveChargeDensityType
struct EffectiveChargeDensity{T, N, S, AT} <: AbstractArray{T, N}

Effective charge density needed to calculate the ElectricPotential. The effective charge density is the charge density (in C/m³) times the volume of the voxel of the respective grid point (in m³). Thus, the unit of the effective charge density is Coulomb (C).

Parametric types

  • T: Element type of data.
  • N: Dimension of the grid and data array.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • data::Array{T, N}: Array containing the values of the effective charge density at the discrete points of the grid.
  • grid::Grid{T, N, S, AT}: Grid defining the discrete points at which the electric potential is determined.
SolidStateDetectors.ElectricFieldType
struct ElectricField{T, N, S, AT} <: AbstractArray{T, N}

Electric field of the simulation in units of volt per meter (V/m).

Parametric types

  • T: Element type of grid.axes.
  • N: Dimension of the grid and data array.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • data::Array{<:StaticArray{Tuple{N}, T}, N}: Array containing the field vectors of the electric field at the discrete points of the grid.
  • grid::Grid{T, N, S, AT}: Grid defining the discrete points for which the electric field is determined.
SolidStateDetectors.ElectricFieldChargeDriftModelType
struct ElectricFieldChargeDriftModel{T <: SSDFloat} <: AbstractChargeDriftModel{T}

Charge drift model in which the electrons and holes drift along the electric field with a mobility of ± 1m²/Vs.

This model is the default when no charge drift model is defined in the configuration file.

SolidStateDetectors.ElectricPotentialType
struct ElectricPotential{T, N, S, AT} <: AbstractArray{T, N}

Electric potential of the simulation in units of volt (V).

Parametric types

  • T: Element type of data.
  • N: Dimension of the grid and data array.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • data::Array{T, N}: Array containing the values of the electric potential at the discrete points of the grid.
  • grid::Grid{T, N, S, AT}: Grid defining the discrete points for which the electric potential is determined.
SolidStateDetectors.EventType
mutable struct Event{T <: SSDFloat}

Collection struct for individual events. This (mutable) struct is meant to be used to look at individual events, not to process a huge amount of events.

Fields

  • locations::VectorOfArrays{CartesianPoint{T}}: Vector of the positions of all hits of the event.
  • energies::VectorOfArrays{T}: Vector of energies corresponding to the hits of the event.
  • drift_paths::Union{Vector{EHDriftPath{T}}, Missing}: Calculated drift paths of each hit position.
  • waveforms::Union{Vector{<:Any}, Missing}: Generated signals (waveforms) of the event.
SolidStateDetectors.GeometricalAzimutalAxisWeightsType
struct GeometricalAzimutalAxisWeights{T} <: AbstractGeometricalAxisWeights{T}

GeometricalAzimutalAxisWeights stores certain precalculated (in its construction) fixed values of the azimulal (or polar) axis, ax, of an cylindrical grid on which the field calculation is performed. These precalculated values are further used in the field calculation to calculate the six weights for the neighboring grid points in the SOR.

GeometricalAzimutalAxisWeights has only one field: weights::Array{T}, 2, which is of size (4, length(ax.ticks)). Each row holding the 4 precalculated values for one axis tick.

Axis ticks, t = ax.ticks, and midpoints, mp = midpoints(get_extended_ticks(ax)),
(the middle points between to axis ticks) inbetween:
... t[i-1] --- mp[i] --- t[i] --- mp[i+1] --- t[i+1] ...
are required for the understanding of the precalculated values.

The columns store the following quantities:

  • weights[1, i]: (mp[i+1] - t[i]) / weights[3, i]
  • weights[2, i]: (t[i] - mp[i]) / weights[3, i]
  • weights[3, i]: mp[i+1] - mp[i]
  • weights[4, i]: inv(t[i] - t[i-1])
  • weights[5, i]: normalize(wφ[3, 1:2:end], 1)
  • weights[6, i]: normalize(wφ[3, 2:2:end], 1)

Developer note: This actually seems to be the same as GeometricalCartesianAxisWeights. Thus, we actually could remove that. Or change this into an alias.

SolidStateDetectors.GeometricalCartesianAxisWeightsType
struct GeometricalCartesianAxisWeights{T} <: AbstractGeometricalAxisWeights{T}

GeometricalCartesianAxisWeights stores certain precalculated (in its construction) fixed values of one of the cartesian (linear) axes, ax, of the grid on which the field calculation is performed. These precalculated values are further used in the field calculation to calculate the six weights for the neighboring grid points in the SOR.

GeometricalCartesianAxisWeights has only one field: weights::Array{T}, 2, which is of size (4, length(ax.ticks)). Each row holding the 4 precalculated values for one axis tick.

Axis ticks, t = ax.ticks, and midpoints, mp = midpoints(get_extended_ticks(ax)),
(the middle points between to axis ticks) inbetween:
... t[i-1] --- mp[i] --- t[i] --- mp[i+1] --- t[i+1] ...
are required for the understanding of the precalculated values.

The columns store the following quantities:

  • weights[1, i]: (mp[i+1] - t[i]) / weights[3, i]
  • weights[2, i]: (t[i] - mp[i]) / weights[3, i]
  • weights[3, i]: mp[i+1] - mp[i]
  • weights[4, i]: inv(t[i] - t[i-1])
SolidStateDetectors.GeometricalRadialAxisWeightsType
struct GeometricalRadialAxisWeights{T} <: AbstractGeometricalAxisWeights{T}

GeometricalRadialAxisWeights stores certain precalculated (in its construction) fixed values of the radial axis, ax, of an cylindrical grid on which the field calculation is performed. These precalculated values are further used in the field calculation to calculate the six weights for the neighboring grid points in the SOR.

GeometricalRadialAxisWeights has only one field: weights::Array{T}, 2, which is of size (6, length(ax.ticks)). Each row holding the 6 precalculated values for one axis tick.

Axis ticks, t = ax.ticks, and midpoints, mp = midpoints(get_extended_ticks(ax)),
(the middle points between to axis ticks) inbetween:
... t[i-1] --- mp[i] --- t[i] --- mp[i+1] --- t[i+1] ...
are required for the understanding of the precalculated values.

The columns store the following quantities:

  • weights[1, i]: (mp[i+1] - t[i]) / weights[3, i]
  • weights[2, i]: (t[i] - mp[i]) / weights[3, i]
  • weights[3, i]: (mp[i+1] - mp[i]) / t[i]
  • weights[4, i]: mp[i+1] / (t[i+1] - t[i])
  • weights[5, i]: mp[i] / (t[i] - t[i-1])
    # Developer note: maybe this one could be removed as it basically stores the same as weights[4, i-1]
  • weights[6, i]: (mp[i+1]^2 - mp[i]^2)/2

Developer note: Some of the weights for i == 1 are set manually as r = 0 requires some special treatment.

SolidStateDetectors.GridType
struct Grid{T, N, S <: AbstractCoordinateSystem, AT} <: AbstractGrid{T, N}

Collection of N axes that define the dimensions of the grid needed to calculate ElectricPotential, ElectricField or WeightingPotential.

Parametric types

  • T: Tick type (element type) of the axes.
  • N: Dimension of the grid.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • axes::AT: Tuple of length N containing DiscreteAxis for each dimension of the grid.

See also DiscreteAxis.

SolidStateDetectors.GridMethod
Grid(sim::Simulation{T, Cartesian}; kwargs...)
Grid(sim::Simulation{T, Cylindrical}; kwargs...)

Initializes a Grid based on the objects defined in a Simulation.

The important points of all objects are sampled and added to the ticks of the grid. The grid initialization can be tuned using a set of keyword arguments listed below.

Arguments

  • sim::Simulation{T, S}: Simulation for which the grid will be defined.

Keywords

  • max_tick_distance = missing: Maximum distance between neighbouring ticks of the grid. Additional grid ticks will be added if two neighbouring ticks are too far apart. max_tick_distance can either be a Quantity, e.g. 1u"mm", or a Tuple of Quantity, e.g. (1u"mm", 15u"°", 3u"mm"), to set it for each axis of the Grid separately. Note that a CartesianGrid3D requires a Tuple{LengthQuantity, LengthQuantity, LengthQuantity} while a CylindricalGrid requires a Tuple{LengthQuantity, AngleQuantity, LengthQuantity}. If max_tick_distance is missing, one fourth of the axis length is used.
  • max_distance_ratio::Real = 5: If the ratio between a tick and its left and right neighbour is greater than max_distance_ratio, additional ticks are added between the ticks that are further apart. This prevents the ticks from being too unevenly spaced.
  • add_ticks_between_important_ticks::Bool = true: If set to true, additional points will be added in between the important points obtained from sampling the objects of the simulation. If some objects are too close together, this will ensure a noticeable gap between them in the calculation of potentials and fields.
  • for_weighting_potential::Bool = false: Grid will be optimized for the calculation of an ElectricPotential if set to true, and of a WeightingPotential if set to false.
SolidStateDetectors.ImpurityScaleType
struct ImpurityScale{T, N, S, AT} <: AbstractArray{T, N}

Impurity scalar field of the simulation. It is kinda an alpha map for the impurity density of semiconductors. It can takes values between 0 and 1:

  • 1: The impurity density has its full value. For grid points in depleted region of the semiconductor.

  • ]0,1[: The impurity density is scaled down but not zero. For grid points at the edge of the depleted region which are partially depleted.

  • 0: The impurity density is set to 0. For grid points in undepleted regions of the semiconductor.

Parametric types

  • T: Element type of data.
  • N: Dimension of the grid and data array.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • data::Array{T, N}: Array containing the values of the impurity scale at the discrete points of the grid.
  • grid::Grid{T, N, S, AT}: Grid defining the discrete points for which the impurity scale is determined.
SolidStateDetectors.LinearChargeDensityType
struct LinearChargeDensity{T <: SSDFloat} <: AbstractChargeDensity{T}

Charge density model which assumes a linear gradient in charge density in each dimension of a Cartesian coordinate system.

Fields

  • offsets::NTuple{3,T}: charge density values at the origin of each dimension.
  • gradients::NTuple{3,T}: linear slopes in x, y and z direction.

Definition in Configuration File

A LinearChargeDensity is defined in the configuration file through the field charge_density (of a passive or surrounding) with name: linear and optional fields x, y and z that can each contain init for initial values at 0 and gradient for gradients in that dimension.

An example definition of a linear charge density looks like this:

charge_density:
  name: cylindrical
  x:  # impurity profile with linear gradient in x
    init: 1.0e-10     # C/m³
    gradient: 1.0e-11 # C/m⁴
SolidStateDetectors.LinearImpurityDensityType
struct LinearImpurityDensity{T <: SSDFloat} <: AbstractImpurityDensity{T}

Impurity density model which assumes a linear gradient in impurity density in each dimension of a Cartesian coordinate system.

Fields

  • offsets::NTuple{3,T}: impurity density values at the origin of each dimension.
  • gradients::NTuple{3,T}: linear slopes in x, y and z direction.

Definition in Configuration File

A LinearImpurityDensity is defined in the configuration file through the field impurity_density (of a passive or surrounding) with name: linear and optional fields x, y and z that can each contain init for initial values at 0 and gradient for gradients in that dimension.

An example definition of a linear impurity density looks like this:

impurity_density:
  name: cylindrical
  x:  # impurity profile with linear gradient in x
    init: 1.0e10     # 1/m³
    gradient: 1.0e11 # 1/m⁴
SolidStateDetectors.NBodyChargeCloudType
struct NBodyChargeCloud{T, N, SH} <: AbstractChargeCloud

Struct which defines a charge cloud consisting of multiple point-like charge carriers, initially distributed around a given center.

Parametric Types

  • T: Precision type.
  • N: Number of shells.
  • SH: Geometry of the shells.

Fields

  • locations::Vector{CartesianPoint{T}}: Positions of the charge carriers that are part of the charge cloud.
  • energies::Vector{T}: Energies of the respective charge carriers, in the same order as locations.
  • shell_structure::SH: Initial geometry of the charge carriers around the center point, relevant for plotting.
SolidStateDetectors.NBodyChargeCloudMethod
NBodyChargeCloud(center::CartesianPoint{T}, energy::T, N::Integer, particle_type::Type{PT} = Gamma; kwargs...)

Returns an NBodyChargeCloud for a given energy deposition at a position that defines the center of the charge cloud, given by a center charge surrounded by shells of approximately N point charges equally distributed on the surface of a sphere. Find the algorithm to create the shells here.

Arguments

  • center::CartesianPoint{T}: Center position of the NBodyChargeCloud.
  • energy::RealQuantity: Deposited energy with units. If no units are given, the value is parsed in units of eV.
  • N::Integer: Approximate number of charges in the NBodyChargeCloud (might vary by around 1%).
  • particle_type: ParticleType of the particle that deposited the energy. Default is Gamma.

Keywords

  • radius::T: Estimate for the radius of the NBodyChargeCloud. Default is determined from particle_type via radius_guess.
  • number_of_shells::Int: Number of shells around the center point. Default is 2.

Example

center = CartesianPoint{T}([0,0,0])
energy = 1460u"keV"
NBodyChargeCloud(center, energy, 200, number_of_shells = 3)
Note

Using values with units for energy requires the package Unitful.jl.

SolidStateDetectors.NBodyChargeCloudMethod
NBodyChargeCloud(center::CartesianPoint{T}, energy::T, particle_type::Type{PT} = Gamma; kwargs...)

Returns an NBodyChargeCloud for a given energy deposition at a position that defines the center of the charge cloud, given by a center charge surrounded by shells consisting of platonic solids.

Arguments

  • center::CartesianPoint{T}: Center position of the NBodyChargeCloud.
  • energy::RealQuantity: Deposited energy with units. If no units are given, the value is parsed in units of eV.
  • particle_type: ParticleType of the particle that deposited the energy. Default is Gamma.

Keywords

  • radius::T: Estimate for the radius of the NBodyChargeCloud. Default is determined from particle_type via radius_guess.
  • number_of_shells::Int: Number of shells around the center point. Default is 2.
  • shell_structure: Geometry with which the charges are distributed in the shells. Default is Dodecahedron.

Example

center = CartesianPoint{T}([0,0,0])
energy = 1460u"keV"
NBodyChargeCloud(center, energy, number_of_shells = 3, shell_structure = SolidStateDetectors.Icosahedron)
Note

Using values with units for energy requires the package Unitful.jl.

SolidStateDetectors.ParticleTypeType
abstract type ParticleType

Type of a particle that deposits energy in a Semiconductor. Currently defined are Alpha, Beta and Gamma.

ParticleType is used to determine the radius of an NBodyChargeCloud, where the default radius for Alpha is 0.1mm and the default radius for Beta and Gamma is 0.5mm.

SolidStateDetectors.PassiveType
struct Passive{T,G,MT,CDM} <: AbstractPassive{T}

Passive object, assigned to a SolidStateDetector.

For the calculation of the ElectricPotential and WeightingPotential, passives can be fixed to a constant potential. They can additionally have a charge density profile that has an influence on the ElectricPotential.

Parametric types

  • T: Precision type.
  • G: Type of geometry.
  • MT: Type of material.
  • CDM: Type of charge_density_model.

Fields

  • name::String: Custom name for the passive, relevant for plotting.
  • id::Int: Unique id that will unambiguously identify the passive.
  • potential::T: Potential (in V) to which the passive will be fixed during the calculation of the electric potential. For floating passives, the potential value is NaN.
  • temperature::T: Temperature (in K) of the passive.
  • material::MT: Material of the passive.
  • charge_density_model::CDM: Charge density model for the points inside the passive.
  • geometry::G: Geometry of the passive, see Constructive Solid Geometry (CSG).

Definition in Configuration File

A Passive is defined through an entry in the passives array of a detector or an entry in the surroundings array in the configuration file. It needs material and geometry and can optionally be given a name, id, potential, temperature and charge_density.

An example definition of passives looks like this:

passives:
  - name: Passivated Surface
    material: HPGe
    charge_density: # ...
    geometry: # ...
  - name: Cryostat
    id: 3
    potential: 0
    temperature: 293K
    material: Al
    geometry: # ...
SolidStateDetectors.PointChargeType
struct PointCharge{T} <: AbstractChargeCloud

Struct which defines a single point-like charge carrier.

Fields

  • locations::SVector{1, CartesianPoint{T}}: Position of the charge carrier, saved as single entry of a Vector.

Example

center = CartesianPoint{T}(0,0,0)
pc = PointCharge(center) # Constructor: creates a PointCharge around (0,0,0)
pc.locations             # Array that only contains the center point
SolidStateDetectors.PointTypeType
const PointType = UInt8

Stores certain information about a grid point via bit-flags.

Right now, there are:

const update_bit      = 0x01
const undepleted_bit  = 0x02
const pn_junction_bit = 0x04
const bulk_bit     = 0x08

Examples

How to get information out of a PointType variable point_type:

  1. point_type & update_bit == 0 -> do not update this point (for fixed points)
  2. point_type & update_bit > 0 -> do update this point
  3. point_type & undepleted_bit > 0 -> this point is undepleted
  4. point_type & pn_junction_bit > 0 -> this point belongs to the solid state detector, meaning that it is in the volume of the n-type or p-type material.
  5. point_type & bulk_bit > 0 -> this point is only surrounded by points marked as pn_junction_bit
SolidStateDetectors.PointTypesType
struct PointTypes{T, N, S, AT} <: AbstractArray{T, N}

Information about the grid points used to calculate the ElectricPotential stored via bit-flags. Data is stored as PointType which is an UInt8.

Parametric types

  • T: Element type of grid.axes.
  • N: Dimension of the grid and data array.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • data::Array{PointType, N}: Array containing the point type values at the discrete points of the grid.
  • grid::Grid{T, N, S, AT}: Grid defining the discrete points for which the point types are determined.

See also PointType.

SolidStateDetectors.PotentialCalculationSetupType
struct PotentialCalculationSetup

PotentialCalculationSetup holds the grid, fields and certain precalculated fixed parameters for the field calculation. This struct will be calculated after each refinement as is depends on the grid.

grid: The 3-dimensional grid, either Cartesian or cylindrical, on which the field will be calculated.

The fields:

  • potential: This is the 4-dimensional array of the extended 3D potential array.

The fourth dimensions comes from the red-black (even-odd) division in order to parallelize the field calculation. Extended means that the grid holds one additional tick at both sides of all axes necessary for boundary handling (reflecting, fixed, ...) at the ends of the grid.

  • point_types: Also a 4-dimensional array. Same structure as the potential-array.
  • volume_weights: Also a 4-dimensional array. Same structure as the potential-array.
  • q_eff_imp: Also a 4-dimensional array. Same structure as the potential-array.
  • imp_scale: Also a 4-dimensional array. Same structure as the potential-array.
  • q_eff_fix: Also a 4-dimensional array. Same structure as the potential-array.
  • ϵ_r: Normal 3-dimensional array! In order to calculate the final 6 weights for each grid point in the SOR,

the eight values of the dielectric permittivity (of each octant) around the grid point needs to be loaded. Here not special division is possible for the red-black (even-odd) division.

Precalculated parameters:

  • geom_weights: The parts of the calculation of the six weights in the SOR can be precalculated. Those are stored here for each axis/dimension.
  • sor_const: Vector holding the SOR constants. In the cartesian case only the first entry is used. As the optimal value for the SOR constant

depends on the grid, the constant is linear increased and the array holds the respective value for each radial axis tick.

  • bias_voltage: maximum_applied_potential - minimum_applied_potential. Used for depletion handling, but might be obsolete by now.
  • maximum_applied_potential: Used for depletion handling, but might be obsolete by now.
  • minimum_applied_potential: Used for depletion handling, but might be obsolete by now.
  • grid_boundary_factors: Used in the application of boundary conditions in the field calculation for decaying (infinite) boundary conditions

to approximate the decay of the potential (depending on the grid).

SolidStateDetectors.SSDIntervalType
struct SSDInterval{T <: SSDFloat, L, R, BL, BR} <: IntervalSets.TypedEndpointsInterval{L,R,T}

Interval containing boundary conditions of left and right boundary as parametric type (BL and BR).

Parametric types

  • T: Precision type.
  • L: Boundary type of the left endpoint.
  • R: Boundary type of the right endpoint.
  • BL: Boundary condition at the left endpoint.
  • BR: Boundary condition at the right endpoint.

The boundary types of an SSDInterval can be L, R ∈ {:closed, :open}.

The boundary conditions of an SSDInterval can be BL, BR ∈ {:periodic, :reflecting, :infinite, :r0, :fixed}.

Fields

  • left::T: Value of the left endpoint.
  • right::T: Value of the right endpoint.
SolidStateDetectors.SemiconductorType
struct Semiconductor{T,G,MT,CDM,IDM} <: AbstractSemiconductor{T}

Semiconductor bulk of a SolidStateDetector.

This is the volume in which electrons and holes will drift during the signal development.

Parametric types

  • T: Precision type.
  • G: Type of geometry.
  • MT: Type of material.
  • CDM: Type of charge_drift_model.
  • IDM: Type of impurity_density_model.

Fields

  • temperature::T: Temperature (in K) of the semiconductor.
  • material::MT: Material of the semiconductor.
  • impurity_density_model::IDM: Impurity density model for the points inside the semiconductor.
  • charge_drift_model::CDM: Model that describes the drift of electrons and holes inside the semiconductor.
  • geometry::G: Geometry of the semiconductor, see Constructive Solid Geometry (CSG).

Definition in Configuration File

A Semiconductor is defined through the semiconductor field of a detector. It needs material and geometry, and can optionally be given a temperature, impurity_density and charge_drift_model.

An example definition of a semiconductor looks like this:

semiconductor:
  material: HPGe
  temperature: 78
  impurity_density: # ...
  charge_drift_model: # ...
  geometry: # ...
SolidStateDetectors.SimulationType
mutable struct Simulation{T <: SSDFloat, CS <: AbstractCoordinateSystem} <: AbstractSimulation{T}

Collection of all parts of a simulation of a SolidStateDetector.

Parametric types

  • T: Precision type.
  • CS: Coordinate system (Cartesian or Cylindrical).

Fields

  • config_dict::Dict: Dictionary (parsed configuration file) which initialized the simulation.
  • input_units::NamedTuple: Units with which the config_dict should be parsed.
  • medium::NamedTuple: Medium of the world.
  • detector::Union{SolidStateDetector{T}, Missing}: The SolidStateDetector of the simulation.
  • world::World{T, 3, CS}: The World of the simulation.
  • q_eff_imp::Union{EffectiveChargeDensity{T}, Missing}: Effective charge resulting from the impurites in the Semiconductor of the detector.
  • imp_scale::Union{ImpurityScale{T}, Missing}: Scale (alpha channel) of the impurity density (for depletion handling).
  • q_eff_fix::Union{EffectiveChargeDensity{T}, Missing}: Fixed charge resulting from fixed space charges in Passive of the detector.
  • ϵ_r::Union{DielectricDistribution{T}, Missing}: The DielectricDistribution of the simulation.
  • point_types::Union{PointTypes{T}, Missing}: The PointTypes of the simulation.
  • electric_potential::Union{ElectricPotential{T}, Missing}: The ElectricPotential of the simulation.
  • weighting_potentials::Vector{Any}: The WeightingPotential for each Contact of the detector in the simulation.
  • electric_field::Union{ElectricField{T}, Missing}: The ElectricField of the simulation.
SolidStateDetectors.SolidStateDetectorType
struct SolidStateDetector{T,SC,CT,PT,VDM} <: AbstractConfig{T}

Struct to describe all parts of a solid state detector, i.e. the Semiconductor, a set of Contact and (optionally) Passive and virtual drift volumes.

The properties of the parts (charge densities, fixed potentials, relative permittivity of the materials) will be used as input to the calculation of ElectricPotential and WeightingPotential in the Simulation.

Parametric types

  • T: Precision type.
  • SC: Type of the semiconductor.
  • CT: Type of the contacts.
  • PT: Type of the passives.
  • VDM: Type of the virtual_drift_volumes.

Fields

  • name::String: Name of the detector.
  • semiconductor::SC: Semiconductor of the detector.
  • contacts::CT: Vector of Contact of the detector.
  • passives::PT: Vector of Passive objects, e.g. holding structures around the detector.
  • virtual_drift_volumes::VDM: Vector of virtual drift volumes in which the drift can be modulated by user-defined methods for modulate_driftvector, e.g. DeadVolume.

See also Semiconductor, Contact, Passive and DeadVolume.

SolidStateDetectors.VelocityParametersType
struct VelocityParameters{T <: SSDFloat}

Values needed to parametrize the longitudinal drift velocity of electrons or hole along a crystal axis as a function of the electric field strength.

Background information

The parameterization for the longitudinal drift velocity, $v_l$, as a function of the electric field strength, $E$, was proposed by D.M. Caughey and R.E. Thomas and later expanded by L. Mihailescu et al.:

\[v_l = \frac{\mu_0 E}{(1 + (E/E_0 )^{\beta})^{1/ \beta}} - \mu_{n} E.\]

with the four parameters, $\mu_0$, $\beta$, $E_0$ and $\mu_n$, which are different for electrons and holes and for the different crystal axes.

Note

The parameter $\mu_n$ accounts for the Gunn effects for electrons and should be 0 for holes.

Fields

  • mu0::T: Parameter $\mu_0$ in the parameterization shown above.
  • beta::T: Parameter $\beta$ in the parameterization shown above.
  • E0::T: Parameter $E_0$ in the parameterization shown above.
  • mun::T: Parameter $\mu_n$ in the parameterization shown above.
SolidStateDetectors.WeightingPotentialType
struct WeightingPotential{T, N, S, AT} <: AbstractArray{T, N}

Weighting potential for a given Contact which is a unitless potential.

Parametric types

  • T: Element type of data.
  • N: Dimension of the grid and data array.
  • S: Coordinate system (Cartesian or Cylindrical).
  • AT: Axes type.

Fields

  • data::Array{T, N}: Array containing the values of the weighting potential at the discrete points of the grid.
  • grid::Grid{T, N, S, AT}: Grid defining the discrete points for which the weighting potential is determined.
SolidStateDetectors.WorldType
struct World{T <: SSDFloat, N, S} <: AbstractWorld{T, N}

Definition of the finite volume on which a Simulation is performed.

Parametric types

  • T: Precision type.
  • N: Dimensions of the world.
  • S: Coordinate system (Cartesian or Cylindrical).

Fields

  • intervals::NTuple{N, SSDInterval{T}}: A set of SSDInterval defining the dimensions of the world.
SolidStateDetectors.add_baseline_and_extend_tailMethod
add_baseline_and_extend_tail(wv::RadiationDetectorSignals.RDWaveform{T,U,TV,UV}, n_baseline_samples::Int, total_waveform_length::Int) where {T,U,TV,UV}

Adds a zero-valued baseline in front of the waveform wv and extends (or cuts off) the waveform at the end with the last value of wv. A waveform of length total_waveform_length is returned.

Arguments

  • wv::RadiationDetectorSignals.RDWaveform{T,U,TV,UV}: A waveform (signal over time).
  • n_baseline_samples::Int: Number of samples added in front of the waveform with values 0.
  • total_waveform_length::Int: Number of samples of the extended waveform which is returned.

Examples

add_baseline_and_extend_tail(wv, 1000, 5000)
Note

This functions assumes that the time steps between the samples of the input waveform wv are the same. Thus, that the input waveform is sampled with a fixed frequency.

SolidStateDetectors.add_fano_noiseFunction
add_fano_noise(E_dep::RealQuantity, E_ionisation::RealQuantity, f_fano::Real)::RealQuantity

Adds Fano noise to an energy deposition E_dep, assuming a detector material ionisation energy E_ionisation and a Fano factor f_fano.

Arguments

  • E_dep::RealQuantity: Energy deposited in a semiconductor material.
  • E_ionisation: Energy needed to create one electron-hole-pair in the semiconductor material.
  • f_fano: Fano factor of the material.

Example

add_fano_noise(100u"keV", 2.95u"eV", 0.129)

Some material properties are stored in SolidStateDetectors.material_properties and can be used here:

material = SolidStateDetectors.material_properties[:HPGe]
add_fano_noise(100u"keV", material.E_ionisation, material.f_fano)
Note

Using values with units for E_dep or E_ionisation requires the package Unitful.jl.

SolidStateDetectors.apply_initial_state!Method
apply_initial_state!(sim::Simulation{T}, ::Type{ElectricPotential}, grid::Grid{T} = Grid(sim);
        not_only_paint_contacts::Bool = true, paint_contacts::Bool = true)::Nothing where {T <: SSDFloat}

Applies the initial state for the calculation of the ElectricPotential. It overwrites sim.electric_potential, sim.q_eff_imp, sim.q_eff_fix, sim.ϵ and sim.point_types with the material properties and fixed potentials defined in sim.detector.

Arguments

  • sim::Simulation{T}: Simulation for which the initial state should be applied.
  • grid::Grid{T}: Grid to apply the initial state on. If no grid is given, a default Grid is determined from sim.

Keywords

  • not_only_paint_contacts::Bool = true: Whether to only use the painting algorithm of the surfaces of Contact without checking if points are actually inside them. Setting it to false should improve the performance but the points inside of Contact are not fixed anymore.
  • paint_contacts::Bool = true: Enable or disable the painting of the surfaces of the Contact onto the grid.

Examples

apply_initial_state!(sim, ElectricPotential, paint_contacts = false)
SolidStateDetectors.apply_initial_state!Method
apply_initial_state!(sim::Simulation{T}, ::Type{WeightingPotential}, contact_id::Int, grid::Grid{T} = Grid(sim))::Nothing

Applies the initial state for the calculation of the WeightingPotential for the [Contact}(@ref) with the id contact_id. It overwrites sim.weighting_potentials[contact_id] with the fixed values on the [Contact}(@ref).

Arguments

  • sim::Simulation{T}: Simulation for which the initial state should be applied.
  • contact_id::Int: The id of the Contact for which the WeightingPotential is to be calculated.
  • grid::Grid{T}: Grid to apply the initial state on. If no grid is given, a default Grid is determined from sim.

Keywords

  • not_only_paint_contacts::Bool = true: Whether to only use the painting algorithm of the surfaces of Contact without checking if points are actually inside them. Setting it to false should improve the performance but the points inside of Contact are not fixed anymore.
  • paint_contacts::Bool = true: Enable or disable the painting of the surfaces of the Contact onto the grid.

Examples

apply_initial_state!(sim, WeightingPotential, 1) # =>  applies initial state for weighting potential of contact with id 1
SolidStateDetectors.calc_new_potential_by_neighbors_3DMethod
function calc_new_potential_by_neighbors_3D(
    volume_weight::T,
    weights::NTuple{6,T},
    neighbor_potentials::NTuple{6,T},
) where {T}

Calculates and returns the new potential value of a grid point given the potential values of the six neighbouring grid points (2 in each dimension) neighbor_potentials, the weights corresponding to those six potentials of the neighboring grid points, the weight due to the volume of the grid point itself (the voxel around it).

For more detailed information see Chapter 5.2.2 "Calculation of the Electric Potential" Eqs. 5.31 to 5.37 in https://mediatum.ub.tum.de/doc/1620954/1620954.pdf.

SolidStateDetectors.calculate_capacitance_matrixMethod
calculate_capacitance_matrix(sim::Simulation{T}; consider_multiplicity::Bool = true) where {T}

Calculates the Maxwell Capacitance N×N-Matrix in units of pF, where N is the number of contacts of sim.detector. The individual elements, $c_{i,j}$, are calculated via calculate_mutual_capacitance(sim::Simulation, (i,j)::Tuple{Int,Int}). The matrix should be symmetric. The difference of C[i,j] and C[j,i] are due to numerical precision in the integration due to the different grids of the two weighting potentials.

Arguments

  • sim::Simulation: Simulation for which the capacitance matrix is calculated.

Keywords

  • consider_multiplicity::Bool = true: Whether symmetries of the system should be taken into account. For example, in case of true coaxial detector center around the origin and calculated on a cartesian grid with the x-axis going from [0, x_max] and the y-axis going from [0, y_max] the multiplicity is 4 and, if consider_multiplicity == true, the returned value is already multiplied by 4.
SolidStateDetectors.calculate_electric_field!Method
calculate_electric_field!(sim::Simulation{T}; n_points_in_φ::Union{Missing, Int} = missing)::Nothing

Calculates the ElectricField from the ElectricPotential stored in sim.electric_potential and stores it in sim.electric_field.

Arguments

  • sim::Simulation{T}: Simulation for which sim.electric_potential has already been calculated.

Keywords

  • n_points_in_φ::Union{Missing, Int}: For a 2D ElectricPotential (cylindrical coordinates and symmetric in φ), sim.electric_potential is extended to n_points_in_φ "layers" in φ in order to calculate a 3D [ElectricField]. If n_points_in_φ is missing, the default value is 36.

Examples

calculate_electric_field!(sim, n_points_in_φ = 32)
Note

This method only works if sim.electric_potential has already been calculated and is not missing.

SolidStateDetectors.calculate_electric_potential!Method
calculate_electric_potential!(sim::Simulation{T}; kwargs...)::Nothing

Calculates the ElectricPotential for a given Simulation sim on an adaptive grid through successive over relaxation and stores it in sim.electric_potential.

There are several keyword arguments which can be used to tune the calculation.

Arguments

Keywords

  • convergence_limit::Real: convergence_limit times the bias voltage sets the convergence limit of the relaxation. The convergence value is the absolute maximum difference of the potential between two iterations of all grid points. Default of convergence_limit is 1e-7 (times bias voltage).
  • refinement_limits: Defines the maximum relative (to applied bias voltage) allowed differences of the potential value of neighbored grid points in each dimension for each refinement.
    • rl::Real -> One refinement with rl equal in all 3 dimensions.
    • rl::Tuple{<:Real,<:Real,<:Real} -> One refinement with rl set individual for each dimension.
    • rl::Vector{<:Real} -> length(l) refinements with rl[i] being the limit for the i-th refinement.
    • rl::Vector{<:Real,<:Real,<:Real}} -> length(rl) refinements with rl[i] being the limits for the i-th refinement.
  • min_tick_distance::Tuple{<:Quantity, <:Quantity, <:Quantity}: Tuple of the minimum allowed distance between two grid ticks for each dimension. It prevents the refinement to make the grid too fine. Default is 1e-5 for linear axes and 1e-5 / (0.25 * r_max) for the polar axis in case of a cylindrical grid.
  • max_tick_distance::Tuple{<:Quantity, <:Quantity, <:Quantity}: Tuple of the maximum allowed distance between two grid ticks for each dimension used in the initialization of the grid. Default is 1/4 of size of the world of the respective dimension.
  • max_distance_ratio::Real: Maximum allowed ratio between the two distances in any dimension to the two neighbouring grid points. If the ratio is too large, additional ticks are generated such that the new ratios are smaller than max_distance_ratio. Default is 5.
  • grid::Grid: Initial grid used to start the simulation. Default is Grid(sim).
  • depletion_handling::Bool: Enables the handling of undepleted regions. Default is false.
  • use_nthreads::Union{Int, Vector{Int}}: If <:Int, use_nthreads defines the maximum number of threads to be used in the computation. Fewer threads might be used depending on the current grid size due to threading overhead. Default is Base.Threads.nthreads(). If <:Vector{Int}, use_nthreads[i] defines the number of threads used for each grid (refinement) stage of the field simulation. The environment variable JULIA_NUM_THREADS must be set appropriately before the Julia session was started (e.g. export JULIA_NUM_THREADS=8 in case of bash).
  • sor_consts::Union{<:Real, NTuple{2, <:Real}}: Two element tuple in case of cylindrical coordinates. First element contains the SOR constant for r = 0. Second contains the constant at the outer most grid point in r. A linear scaling is applied in between. First element should be smaller than the second one and both should be ∈ [1.0, 2.0]. Default is [1.4, 1.85]. In case of Cartesian coordinates, only one value is taken.
  • max_n_iterations::Int: Set the maximum number of iterations which are performed after each grid refinement. Default is 10000. If set to -1 there will be no limit.
  • not_only_paint_contacts::Bool = true: Whether to only use the painting algorithm of the surfaces of Contact without checking if points are actually inside them. Setting it to false should improve the performance but the points inside of Contact are not fixed anymore.
  • paint_contacts::Bool = true: Enable or disable the painting of the surfaces of the Contact onto the grid.
  • verbose::Bool=true: Boolean whether info output is produced or not.

Example

calculate_electric_potential!(sim, refinement_limits = [0.3, 0.1, 0.05], max_distance_ratio = 4, max_n_iterations = 20000)
SolidStateDetectors.calculate_mutual_capacitanceMethod
calculate_mutual_capacitance(sim::Simulation, ij::Tuple{Int, Int}; consider_multiplicity::Bool = true)

Returns the mutual capacitance between the contacts with ID i = ij[1] and j = ij[2]. It is calculated via the weighting potentials of the contacts, $\Phi_i^w(\vec{r})$ and $\Phi_j^w(\vec{r})$:

\[c_{ij} = \epsilon_0 \int_{World} \nabla \Phi_i^w(\vec{r}) ϵ_r(\vec{r}) \nabla \Phi_j^w(\vec{r}) d\vec{r}\]

Note

These are elements of the Mawell Capcitance Matrix. Look up Capacitances for more information.

Note

The electric potential as well as the two weighting potentials of both contacts have to be calculated.

Arguments

  • sim::Simulation: Simulation for which the capacitance matrix is calculated.
  • ij::Tuple{Int,Int}: Tuple of indices of the contacts for which the capacitance should be calculated.

Keywords

  • consider_multiplicity::Bool = true: Whether symmetries of the system should be taken into account. For example, in case of true coaxial detector center around the origin and calculated on a cartesian grid with the x-axis going from [0, x_max] and the y-axis going from [0, y_max] the multiplicity is 4 and, if consider_multiplicity == true, the returned value is already multiplied by 4.
SolidStateDetectors.calculate_stored_energyMethod
calculate_stored_energy(sim::Simulation; consider_multiplicity::Bool = true)

Calculates and returns the energy stored in the ElectricField of a SolidStateDetector in a given Simulation in units of J.

Arguments

  • sim::Simulation{T}: Simulation with sim.detector for which the stored energy is calculated.

Keywords

  • consider_multiplicity::Bool = true: Whether symmetries of the system should be taken into account. For example, in case of true coaxial detector center around the origin and calculated on a cartesian grid with the x-axis going from [0, x_max] and the y-axis going from [0, y_max] the multiplicity is 4 and, if consider_multiplicity == true, the returned value is already multiplied by 4.
SolidStateDetectors.calculate_weighting_potential!Method
calculate_weighting_potential!(sim::Simulation{T}, contact_id::Int; kwargs...)::Nothing

Calculates the WeightingPotential for a Contact with contact_id given Simulation sim on an adaptive grid through successive over relaxation and stores it in sim.weighting_potentials[contact_id].

There are several keyword arguments which can be used to tune the calculation.

Keywords

  • convergence_limit::Real: convergence_limit sets the convergence limit of the relaxation. The convergence value is the absolute maximum difference of the potential between two iterations of all grid points. Default of convergence_limit is 1e-7.
  • refinement_limits: Defines the maximum relative allowed differences of the potential value of neighbored grid points in each dimension for each refinement.
    • rl::Real -> One refinement with rl equal in all 3 dimensions.
    • rl::Tuple{<:Real,<:Real,<:Real} -> One refinement with rl set individual for each dimension.
    • rl::Vector{<:Real} -> length(l) refinements with rl[i] being the limit for the i-th refinement.
    • rl::Vector{<:Real,<:Real,<:Real}} -> length(rl) refinements with rl[i] being the limits for the i-th refinement.
  • min_tick_distance::Tuple{<:Quantity, <:Quantity, <:Quantity}: Tuple of the minimum allowed distance between two grid ticks for each dimension. It prevents the refinement to make the grid too fine. Default is 1e-5 for linear axes and 1e-5 / (0.25 * r_max) for the polar axis in case of a cylindrical grid.
  • max_tick_distance::Tuple{<:Quantity, <:Quantity, <:Quantity}: Tuple of the maximum allowed distance between two grid ticks for each dimension used in the initialization of the grid. Default is 1/4 of size of the world of the respective dimension.
  • max_distance_ratio::Real: Maximum allowed ratio between the two distances in any dimension to the two neighbouring grid points. If the ratio is too large, additional ticks are generated such that the new ratios are smaller than max_distance_ratio. Default is 5.
  • grid::Grid: Initial grid used to start the simulation. Default is Grid(sim).
  • depletion_handling::Bool: Enables the handling of undepleted regions. Default is false. This is an experimental feature: In undepleted regions (determined in calculate_electric_potential!(sim; depletion_handling = true)), the dielectric permittivity of the semiconductor is scaled up to mimic conductive behavior. The scale factor can be tuned via the function scaling_factor_for_permittivity_in_undepleted_region.
  • use_nthreads::Union{Int, Vector{Int}}: If <:Int, use_nthreads defines the maximum number of threads to be used in the computation. Fewer threads might be used depending on the current grid size due to threading overhead. Default is Base.Threads.nthreads(). If <:Vector{Int}, use_nthreads[i] defines the number of threads used for each grid (refinement) stage of the field simulation. The environment variable JULIA_NUM_THREADS must be set appropriately before the Julia session was started (e.g. export JULIA_NUM_THREADS=8 in case of bash).
  • sor_consts::Union{<:Real, NTuple{2, <:Real}}: Two element tuple in case of cylindrical coordinates. First element contains the SOR constant for r = 0. Second contains the constant at the outer most grid point in r. A linear scaling is applied in between. First element should be smaller than the second one and both should be ∈ [1.0, 2.0]. Default is [1.4, 1.85]. In case of Cartesian coordinates, only one value is taken.
  • max_n_iterations::Int: Set the maximum number of iterations which are performed after each grid refinement. Default is 10000. If set to -1 there will be no limit.
  • not_only_paint_contacts::Bool = true: Whether to only use the painting algorithm of the surfaces of Contact without checking if points are actually inside them. Setting it to false should improve the performance but the points inside of Contact are not fixed anymore.
  • paint_contacts::Bool = true: Enable or disable the painting of the surfaces of the Contact onto the grid.
  • verbose::Bool=true: Boolean whether info output is produced or not.

Example

calculate_weighting_potential!(sim, 1, refinement_limits = [0.3, 0.1, 0.05], max_distance_ratio = 4, max_n_iterations = 20000)
SolidStateDetectors.drift_charges!Method
drift_charges!(evt::Event{T}, sim::Simulation{T}; kwargs...)::Nothing where {T <: SSDFloat}

Calculates the electron and hole drift paths for the given Event and Simulation and stores them in evt.drift_paths.

Arguments

  • evt::Event{T}: Event for which the charges should be drifted.
  • sim::Simulation{T}: Simulation which defines the setup in which the charges in evt should drift.

Keywords

  • max_nsteps::Int = 1000: Maximum number of steps in the drift of each hit.
  • Δt::RealQuantity = 5u"ns": Time step used for the drift.
  • diffusion::Bool = false: Activate or deactive diffusion of charge carriers via random walk.
  • self_repulsion::Bool = false: Activate or deactive self-repulsion of charge carriers of the same type.
  • verbose = true: Activate or deactivate additional info output.

Example

drift_charges!(evt, sim, Δt = 1u"ns", verbose = false)
Note

Using values with units for Δt requires the package Unitful.jl.

SolidStateDetectors.estimate_depletion_voltageMethod
estimate_depletion_voltage(sim::Simulation{T},
    Umin::T = minimum(broadcast(c -> c.potential, sim.detector.contacts)),
    Umax::T = maximum(broadcast(c -> c.potential, sim.detector.contacts));
    contact_id::Int = determine_bias_voltage_contact_id(sim.detector),
    tolerance::AbstractFloat = 1e-1,
    verbose::Bool = true) where {T <: AbstractFloat}

Estimates the potential needed to fully deplete the detector in a given Simulation at the Contact with id contact_id by bisection method. The default contact_id is determined automatically via determine_bias_voltage_contact_id(sim.detector). The default searching range of potentials is set by the extrema of contact potentials.

Arguments

  • sim::Simulation{T}: Simulation for which the depletion voltage should be determined.
  • Umin::T: The minimum value of the searching range.
  • Umax::T: The maximum value of the searching range.

Keywords

  • contact_id::Int: The id of the Contact at which the potential is applied.
  • tolerance::AbstractFloat: The acceptable accuracy of results. Default is 1e-1.
  • verbose::Bool = true: Activate or deactivate additional info output. Default is true.

Example

using SolidStateDetectors
sim = Simulation(SSD_examples[:InvertedCoax])
calculate_electric_potential!(sim)
estimate_depletion_voltage(sim)
Note

The accuracy of the result depends on the precision of the initial simulation.

Note

This function performs two 2D or 3D field calculations, depending on sim.
Thus, keep in mind that is might consume some memory.

See also is_depleted.

SolidStateDetectors.get_active_volumeMethod
get_active_volume(point_types::PointTypes{T}) where {T}

Returns an approximation of the active volume of the detector by summing up the cell volumes of all cells marked as depleted.

Arguments

  • point_types::PointTypes{T}: Point types of a Simulation.

Examples

get_active_volume(sim.point_types)
Note

Only φ-symmetries are taken into account.

SolidStateDetectors.get_charge_densityFunction
get_charge_density(cd::AbstractChargeDensity, pt::AbstractCoordinatePoint)

Returns the charge density at a given point, pt, based on the charge density model cd.

Arguments

  • cd::AbstractChargeDensity: The AbstractChargeDensity defining the charge density inside a Passive.
  • pt::AbstractCoordinatePoint: The point at which cd is to be evaluated.
Note

The value returned by get_charge_density is in units of C/m³ (SI units).

SolidStateDetectors.get_electron_and_hole_contributionMethod
get_electron_and_hole_contribution(evt::Event{T}, sim::Simulation{T}, contact_id::Int)

Returns the electron and hole contribution to the waveform of a Contact with a given contact_id of an Event as a NamedTuple with two entries: electron_contribution and hole_contribution.

Arguments

  • evt::Event{T}: Event in which the charges have already been drifted.
  • sim::Simulation{T}: Simulation which defines the setup in which the charges in evt were drifted.
  • contact_id::Int: The id of the Contact for which the waveform should be split into electron and hole contribution.

Example

using Plots
using SolidStateDetector
T = Float32

simulation = Simulation{T}(SSD_examples[:InvertedCoax])
simulate!(simulation)
event = Event([CartesianPoint{T}(0.02,0.01,0.05)])
simulate!(event, simulation)

contact_id = 1
wf = get_electron_and_hole_contribution(evt, sim, contact_id)
Note

The drift paths in evt need to be calculated using drift_charges! before calling this function.

See also plot_electron_and_hole_contribution.

SolidStateDetectors.get_impurity_densityFunction
get_impurity_density(id::AbstractImpurityDensity, pt::AbstractCoordinatePoint)

Returns the impurity density at a given point, pt, based on the impurity density model id.

Arguments

  • id::AbstractImpurityDensity: The AbstractImpurityDensity defining the impurity density inside a Semiconductor.
  • pt::AbstractCoordinatePoint: The point at which id is to be evaluated.
Note

The value returned by get_impurity_density is in units of 1/m³ (SI Units).

SolidStateDetectors.get_signals!Method
get_signals!(evt::Event{T}, sim::Simulation{T}; Δt::RealQuantity = 5u"ns")::Nothing where {T <: SSDFloat}

Generates the signals/waveforms from the drift paths of an Event for each Contact, for which a WeightingPotential is specified in sim.weighting_potentials.

The output is stored in evt.waveforms.

Arguments

  • evt::Event{T}: Event for which the waveforms should be generated.
  • sim::Simulation{T}: Simulation which defines the setup in which the waveforms are generated.

Keywords

  • Δt::RealQuantity = 5u"ns": Time steps with which the drift paths were calculated.

Example

get_signals!(evt, sim, Δt = 1u"ns") # if evt.drift_paths were calculated in time steps of 1ns
Note

This method only works if evt.drift_paths has already been calculated and is not missing.

SolidStateDetectors.get_sor_kernelMethod
function get_sor_kernel(::Type{S}, args...)

where S is either Cartesian or Cylindrical.

Developer notes: Currently (February 2022), there are some limitations to the @kernel macro of the package KernelAbstractions.jl. Especially, regarding usage of dispatch.

Thus, we have to write two kernel functions right now for the Cartesian & Cylindrical case: sor_cyl_gpu! and sor_car_gpu!.

Inside kernel functions, everything is (and has to be) inlined and we can make use of multiple dispatch. So in the end we only have to write one function for the kernel, sor_kernel, which is then inlined inside the two kernel functions.

We can also use most of the CPU functions with the restriction that all types have to be independent on the GPU-indices of the kernel. E.g., making use of i23_is_even_t = Val(iseven(i2 + i3)) and other similar statements is not possible, which are used in the CPU implementation for optimization. On the GPU (currently) those statements have to be calculated and booleans have to be passed. Maybe this will change in the future.

SolidStateDetectors.get_ticks_at_positions_of_large_gradientMethod
get_ticks_at_positions_of_large_gradient(epot::ElectricPotential)

The electric potential is analyzed in order to find and return ticks where the gradient (electric field) is strong (relativ to its maximum).

In the 1D case of a pn-junction, the electric field strength is the largest at the position of the pn-junction. Thus, this function is likely to return ticks which are located close to the pn-junction of a semiconductor.

SolidStateDetectors.handle_depletionMethod
function handle_depletion(
    new_potential::T, 
    imp_scale::T,
    neighbor_potentials::NTuple{6,T}, 
    q_eff_imp::T, 
    volume_weight::T,
)::Tuple{T, PointType} where {T}

This function handles the grid points with volumes which are not fully depleted. The decision depends on the proposal for the new potential value, new_potential_proposal, for the respective grid point in that iteration and the potential values of the neighboring grid points, neighbor_potentials.

If vmin = minimum(neighbor_potentials) < new_potential_proposal < vmax = maximum(neighbor_potentials) => fully depleted => imp_scale = 1

If it undershoots or overshoots, the impurity density in the grid point is scaled down via reducing imp_scale ∈ [0, 1] of this grid point.

This decision is based on the fact that the potential inside a solid-state detector increases monotonically from a p+-contact towards an n+-contact. Thus, there cannot be local extrema.

Note

If a fixed charge impurity is present, e.g. due to a charged passivated surface, this handling is probably not valid anymore as the potential between a p+-contact towards an n+-contact is not required to change monotonically anymore.

SolidStateDetectors.merge_second_order_important_ticksMethod
merge_second_order_important_ticks(imp::Vector{T}, imp_second_order::Vector{T}; min_diff::T = T(1e-6)) where {T}

Merge all elements of the second vector, imp_second_order, into the first vector, imp, if they are not too close (via min_diff) to elements of the first vector. Returns the merged vector sorted.

SolidStateDetectors.parse_config_fileMethod
parse_config_file(filename::AbstractString)::Dict where {T <: SSDFloat}

Reads in a configuration file and returns a parsed dictionary which holds all the information specified in the configuration file.

Find detailed information on configuration files in Configuration Files.

Arguments

  • filename::AbstractString: File name of the configuration file. If the file is not in the same directory, a path to the file is required.
Note

Currently supported formats for the configuration files: - YAML: filename ends with .yaml. - JSON: filename ends with .json. - SigGen: filename ends with .config.

SolidStateDetectors.plot_electron_and_hole_contributionFunction
plot_electron_and_hole_contribution(evt::Event{T}, sim::Simulation{T}, contact_id::Int)

Plots the waveform as well as the electron and hole contribution to the waveform of a Contact with a given contact_id of an Event.

Arguments

  • evt::Event{T}: Event in which the charges have already been drifted.
  • sim::Simulation{T}: Simulation which defines the setup in which the charges in evt were drifted.
  • contact_id::Int: The id of the Contact for which the waveform should be split into electron and hole contribution.

Keywords

  • n_samples::Int: Number of samples with which the waveforms will be plotted. The default is the number of samples of the original waveform.

Example

using Plots
using SolidStateDetector
T = Float32

simulation = Simulation{T}(SSD_examples[:InvertedCoax])
simulate!(simulation)
event = Event([CartesianPoint{T}(0.02,0.01,0.05)])
simulate!(event, simulation)

contact_id = 1
plot_electron_and_hole_contribution(evt, sim, contact_id, n_samples = 300)
Note

The drift paths in evt need to be calculated using drift_charges! before calling this function.

Note

This method requires to load the package Plots.jl.

See also get_electron_and_hole_contribution.

SolidStateDetectors.readsiggenMethod
readsiggen(file_path::String; T::Type)

Reads a SigGen configuration file (ending in .config) in 'file_path' and returns a dictionary of all parameters. Non-existing parameteres are set to 0.

Arguments

  • file_path::String: File path leading to the SigGen configuration file.

Keywords

  • T::Type: Type of the parameters in the output dictionary. Default is Float64.
SolidStateDetectors.refine!Method
refine!(sim::Simulation{T}, ::Type{ElectricPotential}, max_diffs::Tuple, minimum_distances::Tuple, kwargs...)

Takes the current state of sim.electric_potential and refines it with respect to the input arguments max_diffs and minimum_distances by

  1. extending the grid of sim.electric_potential to be a closed grid in all dimensions,
  2. refining the axis of the grid based on max_diffs and minimum_distances: Insert new ticks between two existing ticks such that the potential difference between each tick becomes smaller than max_diff[i] (i -> dimension) but that the distances between the ticks stays larger than minimum_distances[i], and
  3. creating the new data array for the refined grid and fill it by interpolation of the the initial grid.

Arguments

  • sim::Simulation{T}: Simulation for which sim.electric_potential will be refined.
  • max_diffs::Tuple{<:Real,<:Real,<:Real}: Maximum potential difference between two discrete ticks of sim.electric_potential.grid after refinement.
  • minimum_distances::Tuple{<:Real,<:Real,<:Real}: Minimum distance (in SI Units) between two discrete ticks of sim.electric_potential.grid after refinement.

Examples

SolidStateDetectors.refine!(sim, ElectricPotential, max_diffs = (100, 100, 100), minimum_distances = (0.01, 0.02, 0.01))
SolidStateDetectors.refine!Method
refine!(sim::Simulation{T}, ::Type{WeightingPotential}, max_diffs::Tuple{<:Real,<:Real,<:Real}, minimum_distances::Tuple{<:Real,<:Real,<:Real})

Takes the current state of sim.weighting_potentials[contact_id] and refines it with respect to the input arguments max_diffs and minimum_distances by

  1. extending the grid of sim.weighting_potentials[contact_id] to be a closed grid in all dimensions,
  2. refining the axis of the grid based on max_diffs and minimum_distances: Insert new ticks between two existing ticks such that the potential difference between each tick becomes smaller than max_diff[i] (i -> dimension) but that the distances between the ticks stays larger than minimum_distances[i], and
  3. creating the new data array for the refined grid and fill it by interpolation of the the initial grid.

Arguments

  • sim::Simulation{T}: Simulation for which sim.weighting_potentials[contact_id] will be refined.
  • contact_id::Int: The id of the Contact for which the WeightingPotential is refined.
  • max_diffs::Tuple{<:Real,<:Real,<:Real}: Maximum potential difference between two discrete ticks of sim.weighting_potentials[contact_id].grid after refinement.
  • minimum_distances::Tuple{<:Real,<:Real,<:Real}: Minimum distance (in SI Units) between two discrete ticks of sim.weighting_potentials[contact_id].grid after refinement.

Examples

SolidStateDetectors.refine!(sim, WeightingPotential, 1, max_diffs = (0.01, 0.01, 0.01), minimum_distances = (0.01, 0.02, 0.01))
SolidStateDetectors.scaling_factor_for_permittivity_in_undepleted_regionMethod
scaling_factor_for_permittivity_in_undepleted_region(sc::Semiconductor{T})::T where {T}

This function is called in the calculations of weighting potentials of undepleted detectors. The electric permittivity, $ϵ_{r}$, is scaled with this function in areas where the detector is undepleted. A value between [0, +Inf] should be returned. However, Inf should not be returned but instead a very high value should be returned in order to mimic perfect conductivity if that is desired.

Arguments

  • sc::Semiconductor{T}: Semiconductor for which the dielectric permittivity should be scaled up.
Experimental feature!

This feature is under research. The goal is to study the properties / signal response of undepleted detector. This function is indented to be overwritten by the user to study the response.

SolidStateDetectors.siggentodictMethod
siggentodict(config::Dict; units::Dict)

Converts the dictionary containing the parameters from a SigGen configuration file to a dictionary that can be understood by SolidStateDetectors.jl.

Arguments

  • config::Dict: Dictionary containing SigGen parameters (output of readsiggen()`).

Keywords

  • units::Dict: Units used in SigGen configuration file (set to "mm", "deg", "V" and "K"). The dictionary needs the fields "length", "angle", "potential" and "temperature".
SolidStateDetectors.simulate!Method
simulate!( sim::Simulation{T}; kwargs...) where {T, S}

Performs a full chain simulation for a given Simulation by

  1. calculating the ElectricPotential,
  2. calculating the ElectricField,
  3. calculating the WeightingPotential for each Contact.

The output is stored in sim.electric_potential, sim.electric_field and sim.weighting_potentials, respectively.

There are several keyword arguments which can be used to tune the simulation.

Arguments

  • sim::Simulation{T}: Simulation for which the full chain simulation should be performed.

Keywords

  • convergence_limit::Real: convergence_limit times the bias voltage sets the convergence limit of the relaxation. The convergence value is the absolute maximum difference of the potential between two iterations of all grid points. Default of convergence_limit is 1e-7 (times bias voltage).
  • refinement_limits: Defines the maximum relative (to applied bias voltage) allowed differences of the potential value of neighboured grid points in each dimension for each refinement.
    • rl::Real -> One refinement with rl equal in all 3 dimensions.
    • rl::Tuple{<:Real,<:Real,<:Real} -> One refinement with rl set individual for each dimension.
    • rl::Vector{<:Real} -> length(l) refinements with rl[i] being the limit for the i-th refinement.
    • rl::Vector{<:Real,<:Real,<:Real}} -> length(rl) refinements with rl[i] being the limits for the i-th refinement.
  • min_tick_distance::Tuple{<:Quantity, <:Quantity, <:Quantity}: Tuple of the minimum allowed distance between two grid ticks for each dimension. It prevents the refinement to make the grid too fine. Default is 1e-5 for linear axes and 1e-5 / (0.25 * r_max) for the polar axis in case of a cylindrical grid.
  • max_tick_distance::Tuple{<:Quantity, <:Quantity, <:Quantity}: Tuple of the maximum allowed distance between two grid ticks for each dimension used in the initialization of the grid. Default is 1/4 of size of the world of the respective dimension.
  • max_distance_ratio::Real: Maximum allowed ratio between the two distances in any dimension to the two neighbouring grid points. If the ratio is too large, additional ticks are generated such that the new ratios are smaller than max_distance_ratio. Default is 5.
  • depletion_handling::Bool: Enables the handling of undepleted regions. Default is false.
  • use_nthreads::Int: Number of threads to use in the computation. Default is Base.Threads.nthreads(). The environment variable JULIA_NUM_THREADS must be set appropriately before the Julia session was started (e.g. export JULIA_NUM_THREADS=8 in case of bash).
  • sor_consts::Union{<:Real, NTuple{2, <:Real}}: Two element tuple in case of cylindrical coordinates. First element contains the SOR constant for r = 0. Second contains the constant at the outer most grid point in r. A linear scaling is applied in between. First element should be smaller than the second one and both should be ∈ [1.0, 2.0]. Default is [1.4, 1.85]. In case of Cartesian coordinates, only one value is taken.
  • max_n_iterations::Int: Set the maximum number of iterations which are performed after each grid refinement. Default is 10000. If set to -1 there will be no limit.
  • not_only_paint_contacts::Bool = true: Whether to only use the painting algorithm of the surfaces of Contact without checking if points are actually inside them. Setting it to false should improve the performance but the points inside of Contact are not fixed anymore.
  • paint_contacts::Bool = true: Enable or disable the painting of the surfaces of the Contact onto the grid.
  • verbose::Bool=true: Boolean whether info output is produced or not.

See also calculate_electric_potential!, calculate_electric_field! and calculate_weighting_potential!.

Example

simulate!(sim, refinement_limits = [0.3, 0.1, 0.05], max_distance_ratio = 4, max_n_iterations = 20000)
SolidStateDetectors.simulate!Method
simulate!(evt::Event{T}, sim::Simulation{T}; kwargs...)::Nothing where {T <: SSDFloat}

Simulates the waveforms for the Event for a given Simulation by

  1. calculating the drift paths of all energy hits, at evt.locations and
  2. generating the waveforms for each Contact, for which a WeightingPotential is specified in sim.weighting_potentials.

The output is stored in evt.drift_paths and evt.waveforms.

Arguments

  • evt::Event{T}: Event for which the charges should be drifted.
  • sim::Simulation{T}: Simulation which defines the setup in which the charges in evt should drift.

Keywords

  • max_nsteps::Int = 1000: Maximum number of steps in the drift of each hit.
  • Δt::RealQuantity = 5u"ns": Time step used for the drift.
  • diffusion::Bool = false: Activate or deactive diffusion of charge carriers via random walk.
  • self_repulsion::Bool = false: Activate or deactive self-repulsion of charge carriers of the same type.
  • verbose = true: Activate or deactivate additional info output.

Example

simulate!(evt, sim, Δt = 1u"ns", verbose = false)

See also drift_charges! and get_signals!.

SolidStateDetectors.simulate_waveformsMethod
simulate_waveforms( mcevents::TypedTables.Table, sim::Simulation{T}; kwargs...)
simulate_waveforms( mcevents::TypedTables.Table, sim::Simulation{T}, output_dir::AbstractString, output_base_name::AbstractString; kwargs...)

Simulates the waveforms for all events defined in mcevents for a given Simulation by

  1. calculating the drift paths of all energy hits defined in mcevents,
  2. determining the signal (waveforms) for each Contact, for which a WeightingPotential is specified in sim.weighting_potentials.

Arguments

  • mcevents::TypedTables.Table: Table with information about events in the simulated setup.
  • sim::Simulation{T}: Simulation which defines the setup in which the charges in mcevents should drift.

If HDF5.jl is loaded, this function has additional arguments.

Additional Arguments (HDF5)

  • output_dir::AbstractString: Directory where the HDF5 output file is saved.
  • output_base_name::AbstractString: Basename of the HDF5 output file, default is "generated_waveforms".

Keywords

  • max_nsteps::Int = 1000: Maximum number of steps in the drift of each hit.
  • Δt::RealQuantity = 4u"ns": Time step used for the drift.
  • diffusion::Bool = false: Activate or deactive diffusion of charge carriers via random walk.
  • self_repulsion::Bool = false: Activate or deactive self-repulsion of charge carriers of the same type.
  • number_of_carriers::Int = 1: Number of charge carriers to be used in the N-Body simulation of an energy deposition.
  • number_of_shells::Int = 1: Number of shells around the center point of the energy deposition.
  • verbose = false: Activate or deactivate additional info output.
  • chunk_n_physics_events::Int = 1000 (HDF5 only): Number of events that should be saved in a single HDF5 output file.

Examples

simulate_waveforms(mcevents, sim, Δt = 1u"ns", verbose = false)
# => returns the input table `mcevents` with an additional column `waveform` in which the generated waveforms are stored
import HDF5
simulate_waveforms(mcevents, sim, "output_dir", "my_basename", Δt = 1u"ns", verbose = false)
# => simulates the charge drift and saves the output to "output_dir/my_basename_evts_xx.h5"
Note

The drift paths are just calculated temporarily and not returned.

Note

Using values with units for Δt requires the package Unitful.jl.

SolidStateDetectors.ssd_readFunction
ssd_read(filename::AbstractString, ::Type{Simulation})

Reads a Simulation from a HDF5 file with a given filename using LegendHDF5IO.jl.

Arguments

  • filename::AbstractString: Filename of the HDF5 file.

Example

using HDF5 
using LegendHDF5IO
using SolidStateDetectors
sim = ssd_read("example_sim.h5", Simulation)
Note

In order to use this method, the packages HDF5.jl and LegendHDF5IO.jl have to be loaded before loading SolidStateDetectors.jl.

See also ssd_write.

SolidStateDetectors.ssd_writeFunction
ssd_write(filename::AbstractString, sim::Simulation)

Converts a Simulation to a NamedTuple and writes it to a HDF5 file with a given filename using LegendHDF5IO.jl.

Arguments

  • filename::AbstractString: Filename of the HDF5 file.
  • sim::Simulation: Simulation that should be written to the HDF5 file.

Example

using HDF5 
using LegendHDF5IO
using SolidStateDetectors
sim = Simulation(SSD_examples[:InvertedCoax])
simulate!(sim)
ssd_write("example_sim.h5", sim)
Warn

If a file with filename already exists, it will be overwritten by this method.

Note

In order to use this method, the packages HDF5.jl and LegendHDF5IO.jl have to be loaded before loading SolidStateDetectors.jl.

See also ssd_read.

SolidStateDetectors.update!Method
function update!(   
    pcs::PotentialCalculationSetup{T}; 
    ::Nothing, # these two unused arguments are used such that the method
    ::Any;     # is similar to the GPU method for it.
    use_nthreads::Int = Base.Threads.nthreads(), 
    depletion_handling::Val{depletion_handling_enabled} = Val{false}(), 
    is_weighting_potential::Val{_is_weighting_potential} = Val{false}(),
    only2d::Val{only_2d} = Val{false}()

)::Nothing where {T, depletionhandlingenabled, only2d, _isweighting_potential}

This function performs one iteration of the SOR. One iteration consists out of 4 steps:

1) Iterate in parallel over all even points and update their potential. 
2) Apply the boundary conditions at the ends of the grid for all even points. 
3) Iterate in parallel over all odd points and update their potential. 
4) Apply the boundary conditions at the ends of the grid for all odd points.
SolidStateDetectors.update_till_convergence!Method
update_till_convergence!( sim::Simulation{T} ::Type{ElectricPotential}, convergence_limit::Real; kwargs...)::T

Takes the current state of sim.electric_potential and updates it until it has converged.

There are several keyword arguments which can be used to tune the simulation.

Arguments

  • sim::Simulation{T}: Simulation for which sim.electric_potential will be updated.
  • convergence_limit::Real: convergence_limit times the bias voltage sets the convergence limit of the relaxation. The convergence value is the absolute maximum difference of the potential between two iterations of all grid points. Default is 1e-7.

Keywords

  • n_iterations_between_checks::Int: Number of iterations between checks. Default is set to 500.
  • max_n_iterations::Int: Set the maximum number of iterations which are performed after each grid refinement. Default is -1. If set to -1 there will be no limit.
  • depletion_handling::Bool: Enables the handling of undepleted regions. Default is false.
  • use_nthreads::Int: Number of threads to use in the computation. Default is Base.Threads.nthreads(). The environment variable JULIA_NUM_THREADS must be set appropriately before the Julia session was started (e.g. export JULIA_NUM_THREADS=8 in case of bash).
  • not_only_paint_contacts::Bool = true: Whether to only use the painting algorithm of the surfaces of Contact without checking if points are actually inside them. Setting it to false should improve the performance but the points inside of Contact are not fixed anymore.
  • paint_contacts::Bool = true: Enable or disable the painting of the surfaces of the Contact onto the grid.
  • sor_consts::Union{<:Real, NTuple{2, <:Real}}: Two element tuple in case of cylindrical coordinates. First element contains the SOR constant for r = 0. Second contains the constant at the outer most grid point in r. A linear scaling is applied in between. First element should be smaller than the second one and both should be ∈ [1.0, 2.0]. Default is [1.4, 1.85]. In case of Cartesian coordinates, only one value is taken.
  • verbose::Bool=true: Boolean whether info output is produced or not.

Example

SolidStateDetectors.update_till_convergence!(sim, ElectricPotential, 1e-6, depletion_handling = true)
SolidStateDetectors.update_till_convergence!Method
update_till_convergence!( sim::Simulation{T} ::Type{WeightingPotential}, contact_id::Int, convergence_limit::Real; kwargs...)::T

Takes the current state of sim.weighting_potentials[contact_id] and updates it until it has converged.

There are several keyword arguments which can be used to tune the simulation.

Arguments

  • sim::Simulation{T}: Simulation for which sim.weighting_potentials[contact_id] will be updated.
  • contact_id::Int: The id of the Contact for which the WeightingPotential is to be calculated.
  • convergence_limit::Real: convergence_limit times the bias voltage sets the convergence limit of the relaxation. The convergence value is the absolute maximum difference of the potential between two iterations of all grid points. Default is 1e-7.

Keywords

  • n_iterations_between_checks::Int: Number of iterations between checks. Default is set to 500.
  • max_n_iterations::Int: Set the maximum number of iterations which are performed after each grid refinement. Default is -1. If set to -1 there will be no limit.
  • depletion_handling::Bool: Enables the handling of undepleted regions. Default is false. This is an experimental feature: In undepleted regions (determined in calculate_electric_potential!(sim; depletion_handling = true)), the dielectric permittivity of the semiconductor is scaled up to mimic conductive behavior. The scale factor can be tuned via the function scaling_factor_for_permittivity_in_undepleted_region.
  • use_nthreads::Int: Number of threads to use in the computation. Default is Base.Threads.nthreads(). The environment variable JULIA_NUM_THREADS must be set appropriately before the Julia session was started (e.g. export JULIA_NUM_THREADS=8 in case of bash).
  • not_only_paint_contacts::Bool = true: Whether to only use the painting algorithm of the surfaces of Contact without checking if points are actually inside them. Setting it to false should improve the performance but the points inside of Contact are not fixed anymore.
  • paint_contacts::Bool = true: Enable or disable the painting of the surfaces of the Contact onto the grid.
  • sor_consts::Union{<:Real, NTuple{2, <:Real}}: Two element tuple in case of cylindrical coordinates. First element contains the SOR constant for r = 0. Second contains the constant at the outer most grid point in r. A linear scaling is applied in between. First element should be smaller than the second one and both should be ∈ [1.0, 2.0]. Default is [1.4, 1.85]. In case of Cartesian coordinates, only one value is taken.
  • verbose::Bool=true: Boolean whether info output is produced or not.

Example

SolidStateDetectors.update_till_convergence!(sim, WeightingPotential, 1, 1e-6, use_nthreads = 4)