API
Types
SolidStateDetectors.ADLChargeDriftModelSolidStateDetectors.AbstractConfigSolidStateDetectors.AbstractLineSolidStateDetectors.BoxSolidStateDetectors.ConstantChargeDensityModelSolidStateDetectors.ContactSolidStateDetectors.DifferenceSolidStateDetectors.DiscreteAxisSolidStateDetectors.DiscreteAxisSolidStateDetectors.ElectricFieldChargeDriftModelSolidStateDetectors.ElectricPotentialSolidStateDetectors.ElectricPotentialSolidStateDetectors.EventSolidStateDetectors.GeometryUnionSolidStateDetectors.GridSolidStateDetectors.IntersectionSolidStateDetectors.LineSolidStateDetectors.LineSegmentSolidStateDetectors.LinearChargeDensityModelSolidStateDetectors.PointTypeSolidStateDetectors.PointTypesSolidStateDetectors.PointTypesSolidStateDetectors.PotentialSimulationSetupSolidStateDetectors.RaySolidStateDetectors.RectangularCuboidSolidStateDetectors.SimulationSolidStateDetectors.SolidStateDetectorSolidStateDetectors.ZeroChargeDensityModel
Functions
SolidStateDetectors.RBArraySolidStateDetectors.RBExtBy2ArraySolidStateDetectors.RBExtBy2ArraySolidStateDetectors._drift_charge!SolidStateDetectors.add_fano_noiseSolidStateDetectors.apply_initial_state!SolidStateDetectors.apply_initial_state!SolidStateDetectors.calculate_capacitanceSolidStateDetectors.calculate_electric_field!SolidStateDetectors.calculate_electric_potential!SolidStateDetectors.calculate_weighting_potential!SolidStateDetectors.get_active_volumeSolidStateDetectors.get_electron_drift_fieldSolidStateDetectors.get_hole_drift_fieldSolidStateDetectors.get_path_to_example_config_filesSolidStateDetectors.get_rbidx_right_neighbourSolidStateDetectors.innerloops!SolidStateDetectors.innerloops!SolidStateDetectors.nidxSolidStateDetectors.parse_config_fileSolidStateDetectors.point_typeSolidStateDetectors.readsiggenSolidStateDetectors.refine!SolidStateDetectors.refine!SolidStateDetectors.siggentodictSolidStateDetectors.simulate!SolidStateDetectors.update!SolidStateDetectors.update_till_convergence!SolidStateDetectors.update_till_convergence!
Documentation
SolidStateDetectors.ADLChargeDriftModel — TypeADLChargeDriftModel{T <: SSDFloat} <: AbstractChargeDriftModelFields
electrons::CarrierParameters{T}holes::CarrierParameters{T}masses::MassParameters{T}- `phi110::T
gammas::SVector{4, SMatrix{3,3,T}}temperaturemodel::AbstractTemperatureModel{T}
SolidStateDetectors.AbstractConfig — Typeabstract type AbstractConfig{T <: SSDFloat} endSupertype 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.AbstractLine — Typeabstract type AbstractLine{T, N, S} <: AbstractGeometry{T, N} endT: eltype N: Dimension S: Coordinate System: :cartesian or :cylindrical
SolidStateDetectors.Box — Typemutable struct Box{T} <: AbstractGeometry{T, 3}Very simple rectengular box in cartesian coordinates.
SolidStateDetectors.ConstantChargeDensityModel — Typestruct ConstantChargeDensityModel{T <: SSDFloat} <: AbstractChargeDensityModel{T}Returns always a fixed charge density.
SolidStateDetectors.Contact — Typemutable struct Contact{T} <: AbstractContact{T}T: Type of precision.
SolidStateDetectors.Difference — Typestruct Difference{T, N, A, B} <: AbstractSet{T, N}a && !b
SolidStateDetectors.DiscreteAxis — TypeDiscreteAxis{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.DiscreteAxis — MethodDiscreteAxis(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.ElectricFieldChargeDriftModel — TypeElectricFieldChargeDriftModel <: AbstractChargeDriftModelSolidStateDetectors.ElectricPotential — MethodElectricPotential(setup::PotentialSimulationSetup{T, 3, :cartesian} ; kwargs...)::ElectricPotential{T, 3, :cartesian}Extracts the electric potential from setup.
SolidStateDetectors.ElectricPotential — MethodElectricPotential(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.Event — Typemutable 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 hugh amount of events.
SolidStateDetectors.GeometryUnion — Typestruct GeometryUnion{T, N, A, B} <: AbstractSet{T, N}a || b
SolidStateDetectors.Grid — TypeT: tick type
N: N dimensional
S: System (Cartesian, Cylindrical...)SolidStateDetectors.Intersection — Typestruct Intersection{T, N, A, B} <: AbstractSet{T, N}a && b
SolidStateDetectors.Line — Typestruct Line{T, N, S} <: AbstractLine{T, N, S}<––A–––––B––>
SolidStateDetectors.LineSegment — Typestruct LineSegment{T, N, S} <: AbstractLine{T, N, S}
[A----------B]SolidStateDetectors.LinearChargeDensityModel — Typestruct LinearChargeDensityModel{T <: SSDFloat} <: AbstractChargeDensityModel{T}Simple charge density model which assumes a linear gradient in charge density in each dimension. offsets::NTuple{3, T} are the charge densities at 0 and gradients::NTuple{3, T} are the linear slopes in each dimension.
SolidStateDetectors.PointType — Typeconst PointType = UInt8Stores 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:
pt & update_bit == 0-> do not update this point (for fixed points)pt & update_bit > 0-> do update this pointpt & undepleted_bit > 0-> this point is undepletedpt & 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.PointTypes — TypePointTypes{T, N, S} <: AbstractArray{T, N}PointTypes stores the point type of each grid point.
SolidStateDetectors.PointTypes — MethodPointTypes(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.PotentialSimulationSetup — TypePotentialSimulationSetup{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.Ray — Typestruct Ray{T, N, S} <: AbstractLine{T, N, S}
[A----------B---->SolidStateDetectors.RectangularCuboid — Typestruct RectangularCuboid{T} <: AbstractGeometry{T, 3}...
SolidStateDetectors.Simulation — Typemutable struct Simulation{T <: SSDFloat} <: AbstractSimulation{T}Collection of all parts of a Simulation of a Solid State Detector.
SolidStateDetectors.SolidStateDetector — Typemutable struct SolidStateDetector{T <: SSDFloat, CS} <: AbstractConfig{T}CS: Coordinate System: -> :cartesian / :cylindrical
SolidStateDetectors.ZeroChargeDensityModel — Typestruct ZeroChargeDensityModel{T <: SSDFloat} <: AbstractChargeDensityModel{T}Returns always 0.
SolidStateDetectors.RBArray — MethodRBExtBy2Array( et::Type, g::Grid{T, N, :cylindrical} )::Array{et, N + 1} where {T, N}Returns a RedBlack array for the grid g.
SolidStateDetectors.RBExtBy2Array — MethodRBExtBy2Array( 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.RBExtBy2Array — MethodRBExtBy2Array( 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._drift_charge! — Method_drift_charge(...)Before calling this function one should check that startpos is inside det: in(startpos, det
SolidStateDetectors.add_fano_noise — Functionadd_fano_noise(E_dep::RealQuantity, E_ionisation::RealQuantity, f_fano::Real)::RealQuantityAdd 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! — Methodfunction apply_initial_state!(sim::Simulation{T}, ::Type{ElectricPotential}, grid::Grid{T} = Grid(sim.detector))::NothingApplies 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! — Methodfunction apply_initial_state!(sim::Simulation{T}, ::Type{WeightingPotential}, contact_id::Int, grid::Grid{T} = Grid(sim.detector))::NothingApplies 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_capacitance — Methodcalculate_capacitance(sim::Simulation{T})::T where {T <: SSDFloat}Calculates the capacitance of an detector in Farad.
SolidStateDetectors.calculate_electric_field! — Methodcalculate_electric_field!(sim::Simulation{T}, args...; n_points_in_φ::Union{Missing, Int} = missing, kwargs...)::NothingToDo...
SolidStateDetectors.calculate_electric_potential! — Methodcalculate_electric_potential!(sim::Simulation{T}; kwargs...)::NothingCompute 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_limittimes 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 ofconvergence_limitis2e-6(times bias voltage).max_refinements::Int: Number of maximum refinements. Default is2. Set it to0to switch off refinement.refinement_limits::Tuple{<:Real, <:Real, <:Real}: Tuple of refinement limits for each dimension (in case of cylindrical coordinates the order isr,φ,z). A refinement limit (e.g.refinement_limits[1]) times the bias voltage of the detectordetis 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 isGrid(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 isBase.Threads.nthreads(). The environment variableJULIA_NUM_THREADSmust be set appropriately before the Julia session was started (e.g.export JULIA_NUM_THREADS=8in case of bash).sor_consts::Union{<:Real, NTuple{2, <:Real}}: Two element tuple in case of cylindrical coordinates. First element contains the SOR constant forr= 0. Second contains the constant at the outer most grid point inr. 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 is10000. If set to-1there will be no limit.verbose::Bool=true: Boolean whether info output is produced or not.
SolidStateDetectors.calculate_weighting_potential! — Methodcalculate_weighting_potential!(sim::Simulation{T}, contact_id::Int; kwargs...)::NothingCompute 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_limittimes 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 ofconvergence_limitis2e-6(times bias voltage).max_refinements::Int: Number of maximum refinements. Default is2. Set it to0to switch off refinement.refinement_limits::Tuple{<:Real, <:Real, <:Real}: Tuple of refinement limits for each dimension (in case of cylindrical coordinates the order isr,φ,z). A refinement limit (e.g.refinement_limits[1]) times the bias voltage of the detectordetis 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 isGrid(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 isBase.Threads.nthreads(). The environment variableJULIA_NUM_THREADSmust be set appropriately before the Julia session was started (e.g.export JULIA_NUM_THREADS=8in case of bash).sor_consts::Union{<:Real, NTuple{2, <:Real}}: Two element tuple in case of cylindrical coordinates. First element contains the SOR constant forr= 0. Second contains the constant at the outer most grid point inr. 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 is10000. If set to-1there will be no limit.verbose::Bool=true: Boolean whether info output is produced or not.
SolidStateDetectors.get_active_volume — Methodget_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_field — Methodget_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_field — Methodget_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.get_path_to_example_config_files — Methodget_path_to_example_config_files()::StringReturns the path to example config files provided by the package.
SolidStateDetectors.get_rbidx_right_neighbour — Methodget_rbidx_right_neighbour(rbidx::Int, ::Val{true}, ::Val{true})::Intneeds docu...
SolidStateDetectors.innerloops! — Methodinnerloops!( 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},
bulk_is_ptype::Val{_bulk_is_ptype} )::Nothing where {T, even_points, depletion_handling_enabled, _bulk_is_ptype}(Vectorized) inner loop for Cartesian coordinates. This function does all the work in the field calculation.
SolidStateDetectors.innerloops! — Methodinnerloops!( 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},
bulk_is_ptype::Val{_bulk_is_ptype} )::Nothing where {T, even_points, depletion_handling_enabled, _bulk_is_ptype}(Vectorized) inner loop for Cylindrical coordinates. This function does all the work in the field calculation.
SolidStateDetectors.nidx — Methodnidx( rbidx::Int, ::Val{true}, ::Val{true})::Intfirst 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_file — MethodSolidStateDetector{T}(filename::AbstractString)::SolidStateDetector{T} where {T <: SSDFloat}Reads in a config-JSON file and returns an Detector struct which holds all information specified in the config file.
SolidStateDetectors.point_type — MethodFor charge drift...SolidStateDetectors.readsiggen — Methodreadsiggen(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! — Methodfunction 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! — Methodfunction 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.siggentodict — Methodsiggentodict(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! — Methodfunction 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! — Methodupdate!(fssrb::PotentialSimulationSetupRB{T, 3, 4, S}, RBT::DataType)::NothingLoop 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! — Methodfunction update_till_convergence!( sim::Simulation{T} ::Type{ElectricPotential}, convergence_limit::Real; kwargs...)::TTakes the current state of sim.electric_potential and updates it until it has converged.
SolidStateDetectors.update_till_convergence! — Methodfunction update_till_convergence!( sim::Simulation{T} ::Type{WeightingPotential}, contact_id::Int, convergence_limit::Real; kwargs...)::TTakes the current state of sim.weighting_potentials[contact_id] and updates it until it has converged.