Inactive layer simulation: True coaxial detector

Set the physical models for the inactive layer (dead layer) simulation with YAML

  • using drift model for the inactive layer
  • using constant lifetime trapping model
  • using mobility-tied diffusion model (by calling the calculate_mobility function with the drift model)
using Plots
using Unitful
using SolidStateDetectors
T = Float64
Float64

the geometry parameters of the model for following display

mm = T(1/1000)
det_rin = 1*mm
det_z = det_r = 10*mm
z_draw = det_z/2
pn_r = 8.957282*mm # this one was calculated by searching the zero impurity point (displayed in the following section)

sim = Simulation{T}(SSD_examples[:TrueCoaxial])
cfn = SSD_examples[:TrueCoaxial]
print(open(f -> read(f, String), cfn))
plot(sim.detector, xunit = u"mm", yunit = u"mm", zunit = u"mm")
name: TrueCoaxial
units:
  length: mm
  angle: deg
  potential: V
  temperature: K
grid:
  coordinates: cylindrical
  axes:
    r:
      to: 10.01
      boundaries: inf
    phi:
      from: 0
      to: 0
      boundaries:
        left: periodic
        right: periodic
    z:
      from: -0.01
      to: 10.01
      boundaries:
        left: inf
        right: inf
medium: vacuum
detectors:
- semiconductor:
    material: HPGe
    temperature: 90
    impurity_density:
      name: PtypePNjunction
      lithium_annealing_temperature: 623K
      lithium_annealing_time: 18minute
      doped_contact_id: 2
      bulk_impurity_density:
        name: constant
        value: -1e10cm^-3
    charge_drift_model:
      model: InactiveLayerChargeDriftModel
      temperature: 90K
      neutral_impurity_density:
        name: constant
        value: 5.6769e15cm^-3
      # # Uncomment the following lines if you want to use different bulk and surface impurity densities in the charge_drift_model section
      # bulk_impurity_density:
      #   name: constant
      #   value: -1e10cm^-3
      # surface_impurity_density:
      #   name: li_diffusion
      #   lithium_annealing_temperature: 623K
      #   lithium_annealing_time: 18minute
      #   doped_contact_id: 2
    charge_trapping_model:
      model: ConstantLifetime
      parameters:
        τh: 1ms
        τe: 1ms
      model_inactive: ConstantLifetime
      parameters_inactive:
        τh: 1μs
        τe: 1μs
      inactive_layer_geometry:
        tube:
          r:
            from: 8.957282 # set to the depth of the pn junction boundary when lithium diffusion temperature is 623 K, time is 18 minutes
            to: 10.0
          h: 10.0
          origin:
            z: 5.0
    geometry:
      tube:
        r:
          from: 1.0
          to: 10.0
        h: 10.0
        origin:
          z: 5.0
  contacts:
    - name: "P⁺"
      id: 1
      material: HPGe
      potential: -500
      geometry:
        tube:
          r:
            from: 1.0
            to: 1.0
          h: 10.0
          origin:
            z: 5.0
    - name: "N⁺"
      id: 2
      material: HPGe
      potential: 0
      geometry:
        tube:
          r:
            from: 10.0
            to: 10.0
          h: 10.0
          origin:
            z: 5.0

tutorial_det_dl

Display the impurity profile defined

Above position of the PN junction boundary pn_r is calculate by searching the point on the curve where the density is zero.

r_list = (0*mm):(0.01*mm):det_r
imp_list = T[]
for r in r_list
    pt = CylindricalPoint{T}(r, 0, z_draw)
    push!(imp_list, SolidStateDetectors.get_impurity_density(sim.detector.semiconductor.impurity_density_model, pt))
end
imp_list = imp_list ./ 1e6 # in cm^-3
plot(r_list/mm, imp_list, xlabel = "r [mm]", ylabel = "Impurity density [cm\$^{-3}\$]", label = "",
    color = :darkblue, lw = 2, grid = :on, xlims = (0, 10), ylims = (-2e10, 1e11))
vline!([pn_r/mm], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary")

tutorial_imp_dl

Dispaly the mobility curve

using SolidStateDetectors: Electron, Hole
cdm = sim.detector.semiconductor.charge_drift_model
depth_list = 0:(0.01*mm):(det_r-pn_r)
hole_mobility_list = []
electron_mobility_list = []
for depth in depth_list
    r = det_r-depth
    pt = CartesianPoint{T}(r, 0, z_draw)
    push!(hole_mobility_list, SolidStateDetectors.calculate_mobility(cdm, pt, Hole))
    push!(electron_mobility_list, SolidStateDetectors.calculate_mobility(cdm, pt, Electron))
end
plot(depth_list/mm, hole_mobility_list, label = "Hole", lw = 4)
plot!(depth_list/mm, electron_mobility_list, label = "Electron", lw = 4)
plot!(ylabel = "Mobility [cm\$^2\$/V/s]", xlabel = "Depth to surface [mm]",
    legend = :topleft, frame = :box, grid = :on, minorgrid = :on, xticks = 0:0.2:1.2, yticks = 0:0.5:4.5, ylims = [0, 4.5])

tutorial_mob_dl

Calculate the electric field

To simulate the inactive layer, we need to calculate the electric field with very fine grid.

calculate_electric_potential!(sim, max_n_iterations = 10, grid = Grid(sim), verbose = false, depletion_handling = true)
g = sim.electric_potential.grid
ax1, ax2, ax3 = g.axes
bulk_tick_dis = 0.05*mm
dl_tick_dis = 0.01*mm
user_additional_ticks_ax1 = sort(vcat(ax1.interval.left:bulk_tick_dis:pn_r, pn_r:dl_tick_dis:ax1.interval.right))
user_ax1 = typeof(ax1)(ax1.interval, user_additional_ticks_ax1)
user_g = typeof(g)((user_ax1, ax2, ax3))
calculate_electric_potential!(sim, refinement_limits = 0.1, use_nthreads = 8, grid = user_g, depletion_handling = true)
calculate_electric_field!(sim)
calculate_weighting_potential!(sim, 1, use_nthreads = 8, depletion_handling = true)
calculate_weighting_potential!(sim, 2, use_nthreads = 8, depletion_handling = true);
plot(
    begin
        imp = plot(sim.imp_scale, φ = 0, xunit = u"mm", yunit = u"mm", title = "impurity scale")
        vline!([pn_r/mm], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary", legendfontsize = 6)
    end,
    begin
        plot(sim.point_types, φ = 0, xunit = u"mm", yunit = u"mm", title = "point types")
        vline!([pn_r/mm], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary", legendfontsize = 6)
    end,
    begin
        plot(sim.electric_potential, title = "electric potential")
        vline!([pn_r/mm], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary", legendfontsize = 6)
    end,
    begin
        plot(sim.electric_field, title = "electric field", clims = (0, 100*2000))
        vline!([pn_r/mm], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary", legendfontsize = 6)
    end,
    size = (800, 600), layout = (2, 2),
)
┌ Warning: `use_nthreads` was set to `Base.Threads.nthreads() = 1`.
                    The environment variable `JULIA_NUM_THREADS` must be set appropriately before the julia session is started.
@ SolidStateDetectors ~/work/SolidStateDetectors.jl/SolidStateDetectors.jl/src/Simulation/Simulation.jl:887
Simulation: TrueCoaxial
Electric Potential Calculation
Bias voltage: 500.0 V
φ symmetry: Detector is φ-symmetric -> 2D computation.
Precision: Float64
Device: CPU
Max. CPU Threads: 1
Coordinate system: Cylindrical
N Refinements: -> 1
Convergence limit: 1.0e-7  => 5.0e-5 V
Initial grid size: (286, 1, 16)

Grid size: (286, 1, 16) - using 1 threads now
┌ Info: 	In electric field calculation: Keyword `n_points_in_φ` not set.
		Default is `n_points_in_φ = 36`. 2D field will be extended to 36 points in φ.
┌ Warning: `use_nthreads` was set to `Base.Threads.nthreads() = 1`.
                    The environment variable `JULIA_NUM_THREADS` must be set appropriately before the julia session is started.
@ SolidStateDetectors ~/work/SolidStateDetectors.jl/SolidStateDetectors.jl/src/Simulation/Simulation.jl:887
Simulation: TrueCoaxial
Weighting Potential Calculation - ID: 1
φ symmetry: Detector is φ-symmetric -> 2D computation.
Precision: Float64
Device: CPU
Max. CPU Threads: 1
Coordinate system: Cylindrical
N Refinements: -> 3
Convergence limit: 1.0e-7
Initial grid size: (12, 1, 16)

Grid size: (16, 1, 16) - using 1 threads now
Grid size: (24, 1, 16) - using 1 threads now
Grid size: (38, 1, 16) - using 1 threads now
┌ Warning: `use_nthreads` was set to `Base.Threads.nthreads() = 1`.
                    The environment variable `JULIA_NUM_THREADS` must be set appropriately before the julia session is started.
@ SolidStateDetectors ~/work/SolidStateDetectors.jl/SolidStateDetectors.jl/src/Simulation/Simulation.jl:887
Simulation: TrueCoaxial
Weighting Potential Calculation - ID: 2
φ symmetry: Detector is φ-symmetric -> 2D computation.
Precision: Float64
Device: CPU
Max. CPU Threads: 1
Coordinate system: Cylindrical
N Refinements: -> 3
Convergence limit: 1.0e-7
Initial grid size: (12, 1, 16)

Grid size: (16, 1, 16) - using 1 threads now
Grid size: (24, 1, 16) - using 1 threads now
Grid size: (38, 1, 16) - using 1 threads now

tutorial_ef_dl

Pulse shape and CCE simulation

depth_list = (0.1*mm):(0.1*mm):(det_r-pn_r)
pulse_list = []
totTime = 5 # us
totEnergy = 1 # 1 keV --> simulating ~339 carrier pairs
max_nsteps = Int(totTime * 1000)
N = Int(totEnergy*1000÷2.95)
for depth in depth_list
    r = det_r-depth
    energy_depos = fill(T(2.95/1000), N) * u"keV"
    starting_positions = repeat([CartesianPoint{T}(r, deg2rad(0), z_draw)], N)
    evt = Event(starting_positions, energy_depos);
    simulate!(evt, sim, Δt = 1u"ns", max_nsteps = max_nsteps, diffusion = true, end_drift_when_no_field = false, self_repulsion = false)
    charge = ustrip(evt.waveforms[1].signal)
    push!(pulse_list, charge)
end
pulse_plot = plot()
eff_list = []
for (i, depth) in enumerate(depth_list)
    depth = round(depth/mm, digits = 1)
    plot!(pulse_list[i], label = "Depth: $(depth) mm", lw = 2, xlabel = "Time [ns]", ylabel = "Amplitude [e]",
        legend = :topright, grid = :on, minorgrid = :on, frame = :box)
    push!(eff_list, maximum(pulse_list[i])/N)
end
cce_plot = plot(depth_list/mm, eff_list, xlabel = "Depth to surface [mm]", lw = 2, ylabel = "Charge collection efficiency",
    frame = :box, grid = :on, minorgrid = :on, xticks = 0:0.2:1.2, yticks = 0:0.1:1, dpi = 500, color = :black, label = "")
plot(pulse_plot, cce_plot, layout = (1, 2), size = (1000, 400), margin = 5Plots.mm)

tutorial_pulse_cce_dl


This page was generated using Literate.jl.