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;

the geometry parameters of the model for following display

det_rin = 1u"mm"
det_z = det_r = 10u"mm"
z_draw = det_z/2
pn_r = 8.957282u"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

In the configuration file shown above, the position of the PN junction boundary pn_r was calculated as the position where the density is zero.

Display the impurity profile defined

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

tutorial_imp_dl

Display the mobility curve

using SolidStateDetectors: Electron, Hole
cdm = sim.detector.semiconductor.charge_drift_model
depth_list = 0u"mm":0.01u"mm":(det_r-pn_r)
mobility_list = map(depth -> let pt::CartesianPoint{T} = CartesianPoint(det_r - depth, 0, z_draw)
    (µe = SolidStateDetectors.calculate_mobility(cdm, pt, Hole) * 10000u"cm^2/(V*s)",
     µh = SolidStateDetectors.calculate_mobility(cdm, pt, Electron) * 10000u"cm^2/(V*s)")
end, depth_list)
plot(depth_list,  getfield.(mobility_list, :µh), label = "Hole", lw = 4)
plot!(depth_list, getfield.(mobility_list, :µe), label = "Electron", lw = 4)
plot!(xlabel = "Depth to surface / mm", ylabel = "Mobility / cm\$^2\$/Vs", unitformat = :nounit, legend = :topleft, xlims = (0u"mm", det_r-pn_r))

tutorial_mob_dl

Calculate the electric field

To simulate the inactive layer, we need to calculate the electric field on a 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.05u"mm"
dl_tick_dis   = 0.01u"mm"
user_additional_ticks_ax1 = sort(vcat(ax1.interval.left*u"m":bulk_tick_dis:pn_r, pn_r:dl_tick_dis:ax1.interval.right*u"m"))
user_ax1 = typeof(ax1)(ax1.interval, SolidStateDetectors.to_internal_units.(user_additional_ticks_ax1))
user_g = typeof(g)((user_ax1, ax2, ax3))
calculate_electric_potential!(sim, refinement_limits = 0.1, grid = user_g, depletion_handling = true)
calculate_electric_field!(sim)
calculate_weighting_potential!(sim, 1, depletion_handling = true)
calculate_weighting_potential!(sim, 2, depletion_handling = true);
plot(
    begin
        imp = plot(sim.imp_scale, φ = 0, xunit = u"mm", yunit = u"mm", title = "impurity scale")
        vline!([pn_r], 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], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary", legendfontsize = 6)
    end,
    begin
        plot(sim.electric_potential, xunit = u"mm", yunit = u"mm", title = "electric potential")
        vline!([pn_r], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary", legendfontsize = 6)
    end,
    begin
        plot(sim.electric_field, xunit = u"mm", yunit = u"mm", title = "electric field", clims = (0, 100*2000))
        vline!([pn_r], lw = 2, ls = :dash, color = :darkred, label = "PN junction boundary", legendfontsize = 6)
    end,
    size = (800, 600), layout = (2, 2),
)
Simulation: TrueCoaxial
Electric Potential Calculation
Bias voltage: 500.0 V
φ symmetry: Detector is φ-symmetric -> 2D computation.
Precision: Float64
Device: CPU
Max. CPU Threads: 4
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 4 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 φ.
Simulation: TrueCoaxial
Weighting Potential Calculation - ID: 1
φ symmetry: Detector is φ-symmetric -> 2D computation.
Precision: Float64
Device: CPU
Max. CPU Threads: 4
Coordinate system: Cylindrical
N Refinements: -> 3
Convergence limit: 1.0e-7
Initial grid size: (12, 1, 16)

Grid size: (16, 1, 16) - using 4 threads now
Grid size: (24, 1, 16) - using 4 threads now
Grid size: (38, 1, 16) - using 4 threads now
Simulation: TrueCoaxial
Weighting Potential Calculation - ID: 2
φ symmetry: Detector is φ-symmetric -> 2D computation.
Precision: Float64
Device: CPU
Max. CPU Threads: 4
Coordinate system: Cylindrical
N Refinements: -> 3
Convergence limit: 1.0e-7
Initial grid size: (12, 1, 16)

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

tutorial_ef_dl

Pulse shape and CCE simulation

depth_list = 0.1u"mm":0.1u"mm":(det_r-pn_r)
totTime = 5u"µs"
totEnergy = 1u"keV" # --> simulating ~339 carrier pairs
N = Int(totEnergy÷2.95u"eV")

pulse_plot = plot()
eff_list = map(depth -> begin
    r = det_r-depth
    energy_depos = fill(2.95u"eV", N)
    starting_positions = repeat([CartesianPoint(r, 0, z_draw)], N)
    evt = Event(starting_positions, energy_depos)
    simulate!(evt, sim, Δt = 1u"ns", max_nsteps = round(Int, totTime / 1u"ns"), diffusion = true, end_drift_when_no_field = false, self_repulsion = false)
    pulse = evt.waveforms[1]
    plot!(pulse_plot, pulse, label = "Depth: $(round(typeof(depth), depth, digits = 1))", lw = 2, yunit = u"eV/V")
    maximum(pulse.signal)/N
end, depth_list)
plot!(pulse_plot, legend = :topright, xlabel = "Time / ns", ylabel = "Amplitude / e", unitformat = :nounit)
cce_plot = plot(depth_list, eff_list, xlabel = "Depth to surface / mm", ylabel = "Charge collection efficiency", lw = 2, color = :black, label = "", unitformat = :nounit)
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.