Simulation Chain: Inverted Coax Detector
using Plots
using SolidStateDetectors
using Unitful
T = Float32
sim = Simulation{T}(SSD_examples[:InvertedCoax])
plot(sim.detector, xunit = u"mm", yunit = u"mm", zunit = u"mm", size = (700, 700))
By default only the contacts and passives are shown in detector plots. One can also take a look at the semiconductor. The following plots show the semiconductor and slices in the y and z direction.
plot(
plot(sim.detector.semiconductor, st = :samplesurface),
plot(sim.detector.semiconductor, st = :slice, y = 0u"mm", lw = 2),
plot(sim.detector.semiconductor, st = :slice, z = 40u"mm", lw = 2),
layout = (1,3), size = (1500,500), legend = false
)
One can also have a look at how the initial conditions look like on the grid (its starts with a very coarse grid):
apply_initial_state!(sim, ElectricPotential) # optional
plot(
plot(sim.electric_potential), # initial electric potential (boundary conditions)
plot(sim.point_types), # map of different point types: fixed point / inside or outside detector volume / depleted/undepleted
plot(sim.q_eff_imp), # charge density distribution
plot(sim.ϵ_r), # dielectric distribution
layout = (1, 4), size = (1600, 500)
)
Next, calculate the electric potential:
calculate_electric_potential!( sim,
refinement_limits = [0.2, 0.1, 0.05, 0.01])
plot(
plot(sim.electric_potential, φ = 20), # initial electric potential (boundary conditions)
plot(sim.point_types), # map of different point types: fixed point / inside or outside detector volume / depleted/undepleted
plot(sim.q_eff_imp), # charge density distribution
plot(sim.ϵ_r), # dielectric distribution
layout = (1, 4), size = (1600, 500)
)
Simulation: Public Inverted Coax
Electric Potential Calculation
Bias voltage: 3500.0 V
φ symmetry: Detector is φ-symmetric -> 2D computation.
Precision: Float32
Device: CPU
Max. CPU Threads: 1
Coordinate system: Cylindrical
N Refinements: -> 4
Convergence limit: 1.0e-7 => 0.00035 V
Initial grid size: (24, 1, 20)
Grid size: (30, 1, 30) - using 1 threads now:
Grid size: (40, 1, 50) - using 1 threads now:
Grid size: (66, 1, 92) - using 1 threads now:
Grid size: (272, 1, 380) - using 1 threads now:
SolidStateDetectors.jl supports active (i.e. depleted) volume calculation:
is_depleted(sim.point_types)
true
get_active_volume(sim.point_types) # approximation (sum of the volume of cells marked as depleted)
239.12766f0 cm^3
Partially depleted detectors
SolidStateDetectors.jl can also calculate the electric potential of a partially depleted detector:
sim_undep = deepcopy(sim)
sim_undep.detector = SolidStateDetector(sim_undep.detector, contact_id = 2, contact_potential = 500); # V <-- Bias Voltage of Mantle
calculate_electric_potential!( sim_undep,
depletion_handling = true,
convergence_limit = 1e-6,
refinement_limits = [0.2, 0.1, 0.05, 0.01],
verbose = false)
plot(
plot(sim_undep.electric_potential),
plot(sim_undep.point_types),
plot(sim_undep.imp_scale),
layout = (1, 3), size = (1200, 600)
)
is_depleted(sim_undep.point_types)
false
Compare both volumes:
println("Depleted: ", get_active_volume(sim.point_types))
println("Undepleted: ", get_active_volume(sim_undep.point_types));
Depleted: 239.12766f0 cm^3
Undepleted: 156.83362f0 cm^3
The depletion voltage can also be estimated from the simulation.
estimate_depletion_voltage(sim)
1874.7024f0 V
Electric field calculation
Calculate the electric field of the fully depleted detector, given the already calculated electric potential:
calculate_electric_field!(sim, n_points_in_φ = 72)
plot(sim.electric_field, full_det = true, φ = 0.0, size = (700, 700),
xunit = u"mm", yunit = u"mm", zunit = u"V/mm", clims = (0,800).*u"V/mm")
plot_electric_fieldlines!(sim, full_det = true, φ = 0.0)
Simulation of charge drifts
Any charge drift model can be used for the calculation of the electric field. If no model is explicitely given, the ElectricFieldChargeDriftModel
is used. Other configurations are saved in their configuration files and can be found under:
<package_directory>/examples/example_config_files/ADLChargeDriftModel/<config_filename>.yaml.
Set the charge drift model of the simulation:
charge_drift_model = ADLChargeDriftModel()
sim.detector = SolidStateDetector(sim.detector, charge_drift_model)
Now, let's create an "random" multi-site event:
starting_positions = [ CylindricalPoint{T}( 0.020, deg2rad(10), 0.015 ),
CylindricalPoint{T}( 0.015, deg2rad(20), 0.045 ),
CylindricalPoint{T}( 0.022, deg2rad(35), 0.025 ) ]
energy_depos = T[1460, 609, 1000] * u"keV" # are needed later in the signal generation
evt = Event(starting_positions, energy_depos);
time_step = 5u"ns"
drift_charges!(evt, sim, Δt = time_step)
plot(sim.detector, xunit = u"mm", yunit = u"mm", zunit = u"mm", size = (700, 700))
plot!(evt.drift_paths)
Weighting potential calculation
We need weighting potentials to simulate the detector charge signal induced by drifting charges. We'll calculate the weighting potential for the point contact and the outer shell of the detector:
for contact in sim.detector.contacts
calculate_weighting_potential!(sim, contact.id, refinement_limits = [0.2, 0.1, 0.05, 0.01], n_points_in_φ = 2, verbose = false)
end
wp1 = plot(sim.weighting_potentials[1])
plot!(sim.detector, st = :slice, φ = 0)
wp2 = plot(sim.weighting_potentials[2])
plot!(sim.detector, st = :slice, φ = 0)
plot(
wp1, wp2,
size = (900, 700)
)
Detector Capacitance Matrix
After the weighting potentials are calculated, the detector capacitance matrix can be calculated in the Maxwell Capacitance Matrix Notation:
calculate_capacitance_matrix(sim)
2×2 Matrix{Unitful.Quantity{Float32, 𝐈^2 𝐓^4 𝐋^-2 𝐌^-1, Unitful.FreeUnits{(pF,), 𝐈^2 𝐓^4 𝐋^-2 𝐌^-1, nothing}}}:
3.4293 pF -3.4225 pF
-3.42064 pF 5.02959 pF
See Capacitances for more information.
Detector waveform generation
Given an interaction at an arbitrary point in the detector, we can now simulate the charge drift and the time evolution of the charges induced on the contacts (e.g. on the point contact):
simulate!(evt, sim) # drift_charges + signal generation of all channels
p_pc_signal = plot( evt.waveforms[1], lw = 1.5, xlims = (0, 1100), xlabel = "Time", unitformat = :slash,
legend = false, tickfontsize = 12, ylabel = "Charge", guidefontsize = 14)
SolidStateDetectors.jl also allows to separate the observed charge signal into the charge induced by the electrons and the charge induced by the holes.
contact_id = 1
plot_electron_and_hole_contribution(evt, sim, contact_id, xlims = (0, 1100), xlabel = "Time",
legend = :topleft, tickfontsize = 12, ylabel = "Charge", guidefontsize = 14)
This page was generated using Literate.jl.