API

Types

Functions

Documentation

SolidStateDetectors.ADLChargeDriftModelType
ADLChargeDriftModel{T <: SSDFloat, M <: AbstractDriftMaterial, N, TM <: AbstractTemperatureModel{T}} <: AbstractChargeDriftModel{T}

Fields

  • electrons::CarrierParameters{T}
  • holes::CarrierParameters{T}
  • phi110::T
  • γ::SVector{N,SMatrix{3,3,T,9}}
  • parameters::ADLParameters{T}
  • temperaturemodel::TM
SolidStateDetectors.AbstractConfigType
abstract type AbstractConfig{T <: SSDFloat} end

Supertype of all detector/world/object configs.

User defined geometries must be subtype of AbstractConfig{T}.

There are a few functions which must be defined for a user config, e.g. struct UserConfig{T} <: AbstractConfig{T}:

For cylindrical grids:

  • in(pt::cylindrical{T}, config::UserConfig{T})::Bool where {T <: SSDFloat}
  • Grid(config::UserConfig{T})::Grid{T, 3, :cylindrical} where {T <: SSDFloat}
  • get_ρ_and_ϵ(pt::cylindrical{T}, config::UserConfig{T})::Tuple{T, T} where {T <: SSDFloat}
  • set_pointtypes_and_fixed_potentials!(pointtypes::Array{PointType, 3}, potential::Array{T, 3}, grid::Grid{T, 3, :cylindrical}, config::UserConfig{T}; weighting_potential_channel_idx::Union{Missing, Int} = missing)::Nothing where {T <: SSDFloat}

For cartesian grids:

  • in(pt::StaticVector{3, T}, config::UserConfig{T})::Bool
  • Grid(config::UserConfig{T})::Grid{T, 3, :cartesian} where {T <: SSDFloat}
  • get_ρ_and_ϵ(pt::StaticVector{3, T}, config::UserConfig{T})::Tuple{T, T} where {T <: SSDFloat}
  • set_pointtypes_and_fixed_potentials!(pointtypes::Array{PointType, 3}, potential::Array{T, 3}, grid::Grid{T, 3, :cartesian}, config::UserConfig{T}; weighting_potential_channel_idx::Union{Missing, Int} = missing)::Nothing where {T <: SSDFloat}
SolidStateDetectors.AbstractLineType
abstract type AbstractLine{T, N, S} <: AbstractGeometry{T, N} end

T: eltype N: Dimension S: Coordinate System: :cartesian or :cylindrical

SolidStateDetectors.BoxType
mutable struct Box{T} <: AbstractGeometry{T, 3}

Very simple rectengular box in cartesian coordinates.

SolidStateDetectors.CylindricalChargeDensityType
struct CylindricalChargeDensity{T <: SSDFloat} <: AbstractChargeDensity{T}

Simple charge density model which assumes a linear gradient in charge density in each spatial dimension of a cylindrical coordinate system. offsets::NTuple{3, T} are the charge densities at 0 and gradients::NTuple{3, T} are the linear slopes in r and z direction.

SolidStateDetectors.DiscreteAxisType
DiscreteAxis{T, BL, BR} <: AbstractAxis{T, BL, BR}
  • T: Type of ticks
  • BL, BR ∈ {:periodic, :reflecting, :infinite, :r0, :fixed}
  • BL: left boundary condition
  • BR: right boundary condition
  • I: IntervalSets.Interval (closed or open boundaries)
SolidStateDetectors.DiscreteAxisMethod
DiscreteAxis(left_endpoint::T, right_endpoint::T, BL::Symbol, BR::Symbol, L::Symbol, R::Symbol, ticks::AbstractVector{T}) where {T}
  • T: Type of ticks
  • BL, BR ∈ {:periodic, :reflecting, :infinite, :r0, :fixed}
  • L, R {:closed, :open}
  • ticks: Ticks of the axis
SolidStateDetectors.ElectricPotentialMethod
ElectricPotential(setup::PotentialSimulationSetup{T, 3, :cartesian} ; kwargs...)::ElectricPotential{T, 3, :cartesian}

Extracts the electric potential from setup.

SolidStateDetectors.ElectricPotentialMethod
ElectricPotential(setup::PotentialSimulationSetup{T, 3, :cylindrical} ; kwargs...)::ElectricPotential{T, 3, :cylindrical}

Extracts the electric potential from setup and extrapolate it to an 2π grid.

For 2D grids (r and z) the user has to set the keyword n_points_in_φ::Int, e.g.: n_points_in_φ = 36.

SolidStateDetectors.EventType
mutable struct Event{T <: SSDFloat}

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

SolidStateDetectors.LinearChargeDensityType
struct LinearChargeDensity{T <: SSDFloat} <: AbstractChargeDensity{T}

Simple charge density model which assumes a linear gradient in charge density in each dimension of a Cartesian coordinate system. offsets::NTuple{3, T} are the charge densities at 0 and gradients::NTuple{3, T} are the linear slopes in each dimension, x, y and z.

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`

How to get information out of a PointType variable pt:

  1. pt & update_bit == 0 -> do not update this point (for fixed points)
  2. pt & update_bit > 0 -> do update this point
  3. pt & undepleted_bit > 0 -> this point is undepleted
  4. pt & pn_junction_bit > 0 -> this point belongs to the solid state detector. So it is in the volume of the n-type or p-type material.
SolidStateDetectors.PointTypesMethod
PointTypes(setup::PotentialSimulationSetup{T, 3, :cylindrical} ; kwargs...)::PointTypes{T, 3, :cylindrical}

Extracts the electric potential from setup and extrapolate it to an 2π grid.

For 2D grids (r and z) the user has to set the keyword n_points_in_φ::Int, e.g.: n_points_in_φ = 36.

SolidStateDetectors.PotentialSimulationSetupType
PotentialSimulationSetup{T, N, S} <: AbstractPotentialSimulationSetup{T, N}

Collection struct. It holds the grid, the potential, the point types, the charge density and the dielectric distribution.

SolidStateDetectors.SimulationType
mutable struct Simulation{T <: SSDFloat} <: AbstractSimulation{T}

Collection of all parts of a Simulation of a Solid State Detector.

SolidStateDetectors.RBArrayMethod
RBExtBy2Array( et::Type, g::Grid{T, N, :cylindrical} )::Array{et, N + 1} where {T, N}

Returns a RedBlack array for the grid g.

SolidStateDetectors.RBExtBy2ArrayMethod
RBExtBy2Array( et::Type, g::Grid{T, N, :cylindrical} )::Array{et, N + 1} where {T, N}

Returns a RedBlack array for the grid g. The RedBlack array is extended in its size by 2 in each geometrical dimension.

SolidStateDetectors.RBExtBy2ArrayMethod
RBExtBy2Array( et::Type, g::Grid{T, 3, :cartesian} )::Array{et, 4} where {T}

Returns a RedBlack array for the grid g. The RedBlack array is extended in its size by 2 in each geometrical dimension.

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

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

SolidStateDetectors.apply_initial_state!Method
function apply_initial_state!(sim::Simulation{T}, ::Type{ElectricPotential}, grid::Grid{T} = Grid(sim.detector))::Nothing

Applies the initial state of the electric potential calculation. It overwrites sim.electric_potential, sim.ρ, sim.ρ_fix, sim.ϵ and sim.point_types.

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

Applies the initial state of the weighting potential calculation for the contact with the id contact_id. It overwrites sim.weighting_potentials[contact_id].

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

Compute the electric potential for the given Simulation sim on an adaptive grid through successive over relaxation.

There are serveral <keyword arguments> which can be used to tune the computation:

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 2e-6 (times bias voltage).
  • max_refinements::Int: Number of maximum refinements. Default is 2. Set it to 0 to switch off refinement.
  • refinement_limits::Tuple{<:Real, <:Real, <:Real}: Tuple of refinement limits for each dimension (in case of cylindrical coordinates the order is r, φ, z). A refinement limit (e.g. refinement_limits[1]) times the bias voltage of the detector det is the maximum allowed voltage difference between two neighbouring grid points in the respective dimension. When the difference is larger, new points are created inbetween. Default is [1e-5, 1e-5, 1e-5].
  • init_grid_spacing::Tuple{<:Real, <:Real, <:Real}: Tuple of the initial distances between two grid points for each dimension. For normal coordinates the unit is meter. For angular coordinates, the unit is radiance. It prevents the refinement to make the grid to fine. Default is [0.005, 10.0, 0.005]`.
  • min_grid_spacing::Tuple{<:Real, <:Real, <:Real}: Tuple of the mimimum allowed distance between two grid points for each dimension. For normal coordinates the unit is meter. For angular coordinates, the unit is radiance. It prevents the refinement to make the grid to fine. Default is [1e-6, 1e-6, 1e-6].
  • grid::Grid{T, N, S}: Initial grid used to start the simulation. Default is Grid(detector, init_grid_spacing=init_grid_spacing).
  • 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.
  • verbose::Bool=true: Boolean whether info output is produced or not.
SolidStateDetectors.calculate_weighting_potential!Method
calculate_weighting_potential!(sim::Simulation{T}, contact_id::Int; kwargs...)::Nothing

Compute the weighting potential for the contact with id contact_id for the given Simulation sim on an adaptive grid through successive over relaxation.

There are serveral <keyword arguments> which can be used to tune the computation:

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 2e-6 (times bias voltage).
  • max_refinements::Int: Number of maximum refinements. Default is 2. Set it to 0 to switch off refinement.
  • refinement_limits::Tuple{<:Real, <:Real, <:Real}: Tuple of refinement limits for each dimension (in case of cylindrical coordinates the order is r, φ, z). A refinement limit (e.g. refinement_limits[1]) times the bias voltage of the detector det is the maximum allowed voltage difference between two neighbouring grid points in the respective dimension. When the difference is larger, new points are created inbetween. Default is [1e-5, 1e-5, 1e-5].
  • init_grid_spacing::Tuple{<:Real, <:Real, <:Real}: Tuple of the initial distances between two grid points for each dimension. For normal coordinates the unit is meter. For angular coordinates, the unit is radiance. It prevents the refinement to make the grid to fine. Default is [0.005, 10.0, 0.005]`.
  • min_grid_spacing::Tuple{<:Real, <:Real, <:Real}: Tuple of the mimimum allowed distance between two grid points for each dimension. For normal coordinates the unit is meter. For angular coordinates, the unit is radiance. It prevents the refinement to make the grid to fine. Default is [1e-6, 1e-6, 1e-6].
  • grid::Grid{T, N, S}: Initial grid used to start the simulation. Default is Grid(detector, init_grid_spacing=init_grid_spacing).
  • 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.
  • verbose::Bool=true: Boolean whether info output is produced or not.
SolidStateDetectors.get_active_volumeMethod
get_active_volume(pts::PointTypes{T}) where {T}

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

SolidStateDetectors.get_electron_drift_fieldMethod
get_electron_drift_field(ef::Array{SVector{3, T},3}, chargedriftmodel::AbstractChargeDriftModel)::Array{SVector{3,T},3} where {T <: SSDFloat}

Applies the charge drift model onto the electric field vectors. The field vectors have to be in cartesian coordinates.

SolidStateDetectors.get_hole_drift_fieldMethod
get_hole_drift_field(ef::Array{SVector{3, T},3}, chargedriftmodel::AbstractChargeDriftModel)::Array{SVector{3,T},3} where {T <: SSDFloat}

Applies the charge drift model onto the hole field vectors. The field vectors have to be in cartesian coordinates.

SolidStateDetectors.innerloops!Method
innerloops!(  iz::Int, rb_tar_idx::Int, rb_src_idx::Int, gw_x::Array{T, 2}, gw_y::Array{T, 2}, gw_z::Array{T, 2}, fssrb::PotentialSimulationSetupRB{T, 3, 4, :cartesian},
                            update_even_points::Val{even_points},
                            depletion_handling::Val{depletion_handling_enabled},
                            )::Nothing where {T, even_points, depletion_handling_enabled}

(Vectorized) inner loop for Cartesian coordinates. This function does all the work in the field calculation.

SolidStateDetectors.innerloops!Method
innerloops!(  ir::Int, rb_tar_idx::Int, rb_src_idx::Int, gw_r::Array{T, 2}, gw_φ::Array{T, 2}, gw_z::Array{T, 2}, fssrb::PotentialSimulationSetupRB{T, 3, 4, :cylindrical},
                            update_even_points::Val{even_points},
                            depletion_handling::Val{depletion_handling_enabled},
                        )::Nothing where {T, even_points, depletion_handling_enabled}

(Vectorized) inner loop for Cylindrical coordinates. This function does all the work in the field calculation.

SolidStateDetectors.nidxMethod
nidx( rbidx::Int, ::Val{true}, ::Val{true})::Int

first type argument: type of the orgal point (for even points -> Val{true}(), else Val{false}()) second type argument: is sum of other point indices even or odd -> (if sum is even -> Val{true}(), else Val{false}())

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

Reads in a config file and returns an Detector struct which holds all information specified in the config file. Currently supported formats for the config file: .json, .yaml

SolidStateDetectors.readsiggenMethod
readsiggen(file_path::String[, T::Type=Float64])

Read the '*.config' file in 'file_path' for SigGen and returns a dictionary of all parameters. Non-existing parameteres are set to 0. ...

Arguments

  • file_path::String: file path for the SigGen config file.
  • T::Type=Float64: type of the parameters in the output dictionary.

...

SolidStateDetectors.refine!Method
function refine!(sim::Simulation{T}, ::Type{ElectricPotential}, max_diffs::Tuple{<:Real,<:Real,<:Real}, minimum_distances::Tuple{<:Real,<:Real,<:Real})

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

SolidStateDetectors.refine!Method
function 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.

SolidStateDetectors.siggentodictMethod
siggentodict(config::Dict[, units::Dict, detector_name::String])

Converts the dictionary containing the parameters from a SigGen config file to a SSD config dictionary. This dictionary can be saved as a JSON file using the JSON package and 'JSON.print(file, config, 4)'. The 'detector_name' is set to "Public Inverted Coax" by default to inherit the colour scheme. ...

Arguments

  • config::Dict: dictionary containing SigGen parameters (output of readsiggen()).
  • units::Dict: units used in SigGen file (set to 'mm', 'deg', 'V' and 'K').
  • detector_name::String: name of the detector.

...

SolidStateDetectors.simulate!Method
function simulate!(sim::Simulation{T};  max_refinements::Int = 1, verbose::Bool = false,
                                    depletion_handling::Bool = false, convergence_limit::Real = 1e-5 ) where {T <: SSDFloat}

ToDo...

SolidStateDetectors.update!Method
update!(fssrb::PotentialSimulationSetupRB{T, 3, 4, S}, RBT::DataType)::Nothing

Loop over even grid points. A point is even if the sum of its cartesian indicies (of the not extended grid) is even. Even points get the red black index (rbi) = 2. ( -> rbpotential[ inds..., rbi ]).

SolidStateDetectors.update_till_convergence!Method
function 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.

SolidStateDetectors.update_till_convergence!Method
function 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.