Skip to content

Commit

Permalink
Pass charge trapping model to add_signal!
Browse files Browse the repository at this point in the history
  • Loading branch information
fhagemann committed Sep 22, 2024
1 parent 6a2bed2 commit 132e6c5
Show file tree
Hide file tree
Showing 4 changed files with 20 additions and 14 deletions.
5 changes: 3 additions & 2 deletions src/Event/Event.jl
Original file line number Diff line number Diff line change
Expand Up @@ -314,16 +314,17 @@ function get_electron_and_hole_contribution(evt::Event{T}, sim::Simulation{T, S}
signal_e::Vector{T} = zeros(T, length(maximum(map(p -> p.timestamps_e, evt.drift_paths))))
signal_h::Vector{T} = zeros(T, length(maximum(map(p -> p.timestamps_h, evt.drift_paths))))

ctm = sim.detector.semiconductor.charge_trapping_model
for i in eachindex(evt.drift_paths)
energy = flatview(evt.energies)[i]

dp_e::Vector{CartesianPoint{T}} = evt.drift_paths[i].e_path
dp_e_t::Vector{T} = evt.drift_paths[i].timestamps_e
add_signal!(signal_e, dp_e_t, dp_e, dp_e_t, T(-1)*energy, wp, S)
add_signal!(signal_e, dp_e_t, dp_e, dp_e_t, -energy, wp, S, ctm)

dp_h::Vector{CartesianPoint{T}} = evt.drift_paths[i].h_path
dp_h_t::Vector{T} = evt.drift_paths[i].timestamps_h
add_signal!(signal_h, dp_h_t, dp_h, dp_h_t, energy, wp, S)
add_signal!(signal_h, dp_h_t, dp_h, dp_h_t, energy, wp, S, ctm)
end
unitless_energy_to_charge = _convert_internal_energy_to_external_charge(sim.detector.semiconductor.material)
return (electron_contribution = RDWaveform(range(zero(T) * u"ns", step = dt * u"ns", length = length(signal_e)), signal_e * unitless_energy_to_charge),
Expand Down
8 changes: 5 additions & 3 deletions src/MCEventsProcessing/MCEventsProcessing.jl
Original file line number Diff line number Diff line change
Expand Up @@ -64,6 +64,7 @@ function simulate_waveforms( mcevents::TypedTables.Table, sim::Simulation{T};
for contact in contacts if !ismissing(sim.weighting_potentials[contact.id]) push!(contact_ids, contact.id) end end
wpots_interpolated = [ interpolated_scalarfield(sim.weighting_potentials[id]) for id in contact_ids ];
electric_field = interpolated_vectorfield(sim.electric_field)
ctm = sim.detector.semiconductor.charge_trapping_model

unitless_energy_to_charge = _convert_internal_energy_to_external_charge(sim.detector.semiconductor.material)
@info "Detector has $(n_contacts) contact"*(n_contacts != 1 ? "s" : "")
Expand All @@ -78,7 +79,7 @@ function simulate_waveforms( mcevents::TypedTables.Table, sim::Simulation{T};
@info "Generating waveforms..."
waveforms = map(
wpot -> map(
x -> _generate_waveform(x.dps, to_internal_units.(x.edeps), Δt, Δtime, wpot, S, unitless_energy_to_charge),
x -> _generate_waveform(x.dps, to_internal_units.(x.edeps), Δt, Δtime, wpot, S, unitless_energy_to_charge, ctm),
TypedTables.Table(dps = drift_paths, edeps = edeps)
),
wpots_interpolated
Expand Down Expand Up @@ -144,11 +145,12 @@ end


function _generate_waveform( drift_paths::Vector{<:EHDriftPath{T}}, charges::Vector{<:SSDFloat}, Δt::RealQuantity, dt::T,
wpot::Interpolations.Extrapolation{T, 3}, S::CoordinateSystemType, unitless_energy_to_charge) where {T <: SSDFloat}
wpot::Interpolations.Extrapolation{T, 3}, S::CoordinateSystemType, unitless_energy_to_charge,
ctm::AbstractChargeTrappingModel{T} = NoChargeTrappingModel{T}) where {T <: SSDFloat}
timestamps = _common_timestamps( drift_paths, dt )
timestamps_with_units = range(zero(Δt), step = Δt, length = length(timestamps))
signal = zeros(T, length(timestamps))
add_signal!(signal, timestamps, drift_paths, T.(charges), wpot, S)
add_signal!(signal, timestamps, drift_paths, T.(charges), wpot, S, ctm)
RDWaveform( timestamps_with_units, signal * unitless_energy_to_charge)
end

Expand Down
17 changes: 10 additions & 7 deletions src/SignalGeneration/SignalGeneration.jl
Original file line number Diff line number Diff line change
Expand Up @@ -21,10 +21,11 @@ function add_signal!(
pathtimestamps::AbstractVector{T},
charge::T,
wpot::Interpolations.Extrapolation{T, 3},
S::CoordinateSystemType
S::CoordinateSystemType,
ctm::AbstractChargeTrappingModel{T} = NoChargeTrappingModel{T}()
)::Nothing where {T <: SSDFloat}

tmp_signal::Vector{T} = _calculate_signal(NoChargeTrappingModel{T}(), path, pathtimestamps, charge, wpot, S)
tmp_signal::Vector{T} = _calculate_signal(ctm, path, pathtimestamps, charge, wpot, S)
itp = interpolate!( (pathtimestamps,), tmp_signal, Gridded(Linear()))
t_max::T = last(pathtimestamps)
i_max::Int = searchsortedlast(timestamps, t_max)
Expand All @@ -36,15 +37,17 @@ function add_signal!(
end


function add_signal!(signal::AbstractVector{T}, timestamps::AbstractVector{T}, path::EHDriftPath{T}, charge::T, wpot::Interpolations.Extrapolation{T, 3}, S::CoordinateSystemType)::Nothing where {T <: SSDFloat}
add_signal!(signal, timestamps, path.e_path, path.timestamps_e, -charge, wpot, S) # electrons induce negative charge
add_signal!(signal, timestamps, path.h_path, path.timestamps_h, charge, wpot, S)
function add_signal!(signal::AbstractVector{T}, timestamps::AbstractVector{T}, path::EHDriftPath{T}, charge::T,
wpot::Interpolations.Extrapolation{T, 3}, S::CoordinateSystemType, ctm::AbstractChargeTrappingModel{T} = NoChargeTrappingModel{T}())::Nothing where {T <: SSDFloat}
add_signal!(signal, timestamps, path.e_path, path.timestamps_e, -charge, wpot, S, ctm) # electrons induce negative charge
add_signal!(signal, timestamps, path.h_path, path.timestamps_h, charge, wpot, S, ctm)
nothing
end

function add_signal!(signal::AbstractVector{T}, timestamps::AbstractVector{<:RealQuantity}, paths::Vector{<:EHDriftPath{T}}, charges::Vector{T}, wpot::Interpolations.Extrapolation{T, 3}, S::CoordinateSystemType)::Nothing where {T <: SSDFloat}
function add_signal!(signal::AbstractVector{T}, timestamps::AbstractVector{<:RealQuantity}, paths::Vector{<:EHDriftPath{T}}, charges::Vector{T},
wpot::Interpolations.Extrapolation{T, 3}, S::CoordinateSystemType, ctm::AbstractChargeTrappingModel{T} = NoChargeTrappingModel{T}())::Nothing where {T <: SSDFloat}
for ipath in eachindex(paths)
add_signal!(signal, timestamps, paths[ipath], charges[ipath], wpot, S)
add_signal!(signal, timestamps, paths[ipath], charges[ipath], wpot, S, ctm)
end
nothing
end
4 changes: 2 additions & 2 deletions src/Simulation/Simulation.jl
Original file line number Diff line number Diff line change
Expand Up @@ -1210,8 +1210,8 @@ function get_signal(sim::Simulation{T, CS}, drift_paths::Vector{EHDriftPath{T}},
wpot::Interpolations.Extrapolation{T, 3} = interpolated_scalarfield(sim.weighting_potentials[contact_id])
timestamps = _common_timestamps( drift_paths, dt )
signal::Vector{T} = zeros(T, length(timestamps))
add_signal!(signal, timestamps, drift_paths, energy_depositions, wpot, CS)
unitless_energy_to_charge = unitless_energy_to_charge = _convert_internal_energy_to_external_charge(sim.detector.semiconductor.material)
add_signal!(signal, timestamps, drift_paths, energy_depositions, wpot, CS, sim.detector.semiconductor.charge_trapping_model)
unitless_energy_to_charge = _convert_internal_energy_to_external_charge(sim.detector.semiconductor.material)
return RDWaveform( range(zero(T) * unit(Δt), step = T(ustrip(Δt)) * unit(Δt), length = length(signal)), signal * unitless_energy_to_charge)
end

Expand Down

0 comments on commit 132e6c5

Please sign in to comment.