Example 1: Inverted Coax Detector

using Plots
using SolidStateDetectors
using Unitful

T = Float32
sim = Simulation{T}(SSD_examples[:InvertedCoax])

plot(sim.detector, size = (700, 700))

tutorial_det

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,
    Grid(sim, max_tick_distance = 3u"mm")) # 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)
)

tutorial_initial_condition

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
Convergence limit: 1.0e-7 => 0.00035 V
Threads: 1
Coordinate system: Cylindrical
N Refinements: -> 4
Initial grid size: (20, 1, 18)

Grid size: (26, 1, 28)
Grid size: (38, 1, 48)
Grid size: (64, 1, 90)
Grid size: (274, 1, 382)

tutorial_calculated_potential

SolidStateDetectors.jl supports active (i.e. depleted) volume calculation:

get_active_volume(sim.point_types) # approximation (sum of the volume of cells marked as depleted)
239.12335f0 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),
    layout = (1, 2), size = (800, 700)
)
[ Info: Maximum number of iterations reached. (`n_iterations = 51000`)
[ Info: Maximum number of iterations reached. (`n_iterations = 51000`)
[ Info: Maximum number of iterations reached. (`n_iterations = 51000`)
[ Info: Maximum number of iterations reached. (`n_iterations = 51000`)

Checking undepleted regions  20%|████▋                  |  ETA: 0:00:01
Checking undepleted regions 100%|███████████████████████| Time: 0:00:01

Checking undepleted regions  20%|████▋                  |  ETA: 0:00:01
Checking undepleted regions  40%|█████████▎             |  ETA: 0:00:00
Checking undepleted regions  60%|█████████████▊         |  ETA: 0:00:00
Checking undepleted regions  80%|██████████████████▍    |  ETA: 0:00:00
Checking undepleted regions 100%|███████████████████████| Time: 0:00:00

tutorial_calculated_potential_undep

Compare both volumes:

println("Depleted:   ", get_active_volume(sim.point_types))
println("Undepleted: ", get_active_volume(sim_undep.point_types));
Depleted:   239.12335f0 cm^3
Undepleted: 157.1676f0 cm^3

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))
plot_electric_fieldlines!(sim, full_det = true, φ = 0.0)

tutorial_electric_field

Drift field calculation

Given the electric field and a charge drift model, calculate drift fields for electrons and holes. Precalculating the drift fields saves time during charge drift simulation:

Any drift field model can be used for the calculation of the electric field. If no model is explicitely given, the Bruyneel model from the Agata Data Library (ADL) 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)

And apply the charge drift model to the electric field:

calculate_drift_fields!(sim)

Now, let's create an "random" (multiside) 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, size = (700, 700))
plot!(evt.drift_paths)

tutorial_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

plot(
    plot(sim.weighting_potentials[1]),
    plot(sim.weighting_potentials[2]),
    size = (900, 700)
)

tutorial_weighting_potentials

Detector waveform generation

Single-event simulation

Given an interaction at an arbitrary point in the detector, we can now simulate the charge drift and the resulting detector charge signals (e.g. at 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 / ns",
                    legend = false, tickfontsize = 12, ylabel = "Energy / eV", guidefontsize = 14)

tutorial_waveforms


This page was generated using Literate.jl.