API
Modules
Types and constants
RadiationDetectorDSP.AbstractRadFIRFilter
RadiationDetectorDSP.AbstractRadIIRFilter
RadiationDetectorDSP.AbstractRadSigFilter
RadiationDetectorDSP.AbstractRadSigFilterInstance
RadiationDetectorDSP.AbstractSamples
RadiationDetectorDSP.ArrayOfSimilarSamples
RadiationDetectorDSP.BiquadFilter
RadiationDetectorDSP.CPUNormAdaptor
RadiationDetectorDSP.CRFilter
RadiationDetectorDSP.CUSPChargeFilter
RadiationDetectorDSP.ConvolutionFilter
RadiationDetectorDSP.ConvolutionMethod
RadiationDetectorDSP.DNIMethod
RadiationDetectorDSP.DifferentiatorFilter
RadiationDetectorDSP.DirectConvolution
RadiationDetectorDSP.FFTConvolution
RadiationDetectorDSP.FilteringType
RadiationDetectorDSP.FirstOrderIIR
RadiationDetectorDSP.Gauss1DFilter
RadiationDetectorDSP.IntegratorCRFilter
RadiationDetectorDSP.IntegratorFilter
RadiationDetectorDSP.IntegratorModCRFilter
RadiationDetectorDSP.Intersect
RadiationDetectorDSP.InvCRFilter
RadiationDetectorDSP.InvModCRFilter
RadiationDetectorDSP.InvRCFilter
RadiationDetectorDSP.InvSecondOrderCRFilter
RadiationDetectorDSP.LinearFiltering
RadiationDetectorDSP.ModCRFilter
RadiationDetectorDSP.NonlinearFiltering
RadiationDetectorDSP.PolynomialDNI
RadiationDetectorDSP.RCFilter
RadiationDetectorDSP.SamplesOrWaveform
RadiationDetectorDSP.SamplingInfo
RadiationDetectorDSP.SavitzkyGolayFilter
RadiationDetectorDSP.SecondOrderCRFilter
RadiationDetectorDSP.SignalEstimator
RadiationDetectorDSP.SimpleCSAFilter
RadiationDetectorDSP.TrapezoidalChargeFilter
RadiationDetectorDSP.TruncateFilter
RadiationDetectorDSP.ZACChargeFilter
Functions and macros
RadiationDetectorDSP.adapt_memlayout
RadiationDetectorDSP.add_rect_pulse!
RadiationDetectorDSP.bc_rdfilt
RadiationDetectorDSP.bc_rdfilt!
RadiationDetectorDSP.charge_trapflt!
RadiationDetectorDSP.cr_filter
RadiationDetectorDSP.crmod_filter
RadiationDetectorDSP.cusp_charge_filter_coeffs
RadiationDetectorDSP.dfilt!
RadiationDetectorDSP.differentiator_filter
RadiationDetectorDSP.elsmplinfo
RadiationDetectorDSP.flt_input_length
RadiationDetectorDSP.flt_input_smpltype
RadiationDetectorDSP.flt_output_length
RadiationDetectorDSP.flt_output_smpltype
RadiationDetectorDSP.flt_output_time_axis
RadiationDetectorDSP.fltinstance
RadiationDetectorDSP.gaussian_coeffs
RadiationDetectorDSP.gen_rect_pulse
RadiationDetectorDSP.integrator_cr_filter
RadiationDetectorDSP.integrator_crmod_filter
RadiationDetectorDSP.integrator_filter
RadiationDetectorDSP.inv_cr_filter
RadiationDetectorDSP.inv_crmod_filter
RadiationDetectorDSP.inv_rc_filter
RadiationDetectorDSP.multiply_waveform
RadiationDetectorDSP.rc_filter
RadiationDetectorDSP.rdfilt
RadiationDetectorDSP.reverse_waveform
RadiationDetectorDSP.shift_waveform
RadiationDetectorDSP.simple_csa_response_filter
RadiationDetectorDSP.smplinfo
RadiationDetectorDSP.zac_charge_filter_coeffs
Documentation
RadiationDetectorDSP.AbstractRadFIRFilter
— Typeabstract type AbstractRadFIRFilter <: AbstractRadLinearFilter
Abstract type for FIR filters.
RadiationDetectorDSP.AbstractRadIIRFilter
— Typeabstract type AbstractRadIIRFilter <: AbstractRadSigFilter{LinearFiltering}
Abstract type for IIR filters.
RadiationDetectorDSP.AbstractRadSigFilter
— Typeabstract type AbstractRadSigFilter{FT<:FilteringType} <: Function
Abstract type for signal filters.
Filters are callable as (flt::AbstractRadSigFilter)(input)
and come with specialized broadcasting.
Subtypes of AbstractRadSigFilter
must implement
fltinstance(flt::AbstractRadSigFilter, si::SamplingInfo)::[`AbstractRadSigFilterInstance`](@ref)
Invertible filters should also implement
InverseFunctions.inverse(flt::SomeFilter)
Note that while a filter may have an inverse, it may, depending on the filter paramters, be very unstable in the presence of additional noise. Filters with a high-pass characteristic pass high-frequency noise, so their inverses pass such noise as well without amplifying it (substantially). Filters with a low-pass characteristic, on the other hand, attentuate high-frequency noise, so their inverses amplify such noise and are typically not useful to deconvolve signals in practical applications.
RadiationDetectorDSP.AbstractRadSigFilterInstance
— Typeabstract type AbstractRadSigFilterInstance{FT<:FilteringType}
Abstract type for signal filter instances. Filter instances are specilized to a specific length and numerical type of input and output.
Filter instances are callable as (fi::SomeFilterInstance)(input)
and come with specialized broadcasting.
Subtypes of AbstractRadSigFilterInstance
must implement
RadiationDetectorDSP.rdfilt!(output, fi::SomeFilterInstance, input)
RadiationDetectorDSP.flt_output_smpltype(fi::SomeFilterInstance)
RadiationDetectorDSP.flt_input_smpltype(fi::SomeFilterInstance)
RadiationDetectorDSP.flt_output_length(fi::SomeFilterInstance)
RadiationDetectorDSP.flt_input_length(fi::SomeFilterInstance)
RadiationDetectorDSP.flt_output_time_axis(fi::SomeFilterInstance, time::AbstractVector{<:RealQuantity})
Invertible filter instances should implement
InverseFunctions.inverse(fi::SomeFilterInstance)
Default methods are implemented for
RadiationDetectorDSP.rdfilt(fi::AbstractRadSigFilterInstance, x::AbstractSamples)
RadiationDetectorDSP.rdfilt(fi::AbstractRadSigFilterInstance, wf::RDWaveform)
RadiationDetectorDSP.bc_rdfilt(fi::AbstractRadSigFilterInstance, inputs)
The default methods that operate on RadiationDetectorSignals.RDWaveform
s require RadiationDetectorDSP.flt_output_time_axis
.
RadiationDetectorDSP.AbstractSamples
— Typeconst AbstractSamples{T<:RealQuantity} = AbstractVector{T}
A vector of signal samples.
RadiationDetectorDSP.ArrayOfSimilarSamples
— Typeconst ArrayOfSimilarSamples{T<:RealQuantity} = ArrayOfSimilarVectors{T}
An array of similar sample vectors.
RadiationDetectorDSP.BiquadFilter
— Typestruct BiquadFilter{T<:RealQuantity} <: AbstractRadIIRFilter
Constructors:
BiquadFilter(fields...)
Fields:
b_012::Tuple{T, T, T} where T<:Union{Real, Unitful.AbstractQuantity{<:Real}}
: Coefficients b0 to b2a_12::Tuple{T, T} where T<:Union{Real, Unitful.AbstractQuantity{<:Real}}
: Coefficients a1 to a2, a_0 equals one implicitly
RadiationDetectorDSP.CPUNormAdaptor
— TypeRadiationDetectorDSP.CPUNormAdaptor
To be used with Adapt.adapt
.
Adapt.adapt(RadiationDetectorDSP.CPUNormAdaptor, x)
adapts x
to reside on the CPU and tries to ensure that arrays are stored in column-major order.
RadiationDetectorDSP.CRFilter
— Typestruct CRFilter <: AbstractRadIIRFilter
A first-order CR highpass filter.
The inverse filter is InvCRFilter
, this is typically stable even in the presence of additional noise. This is because a CR filter passes high-frequency noise and so it's inverse passes such noise as well without amplifying it.
Constructors:
CRFilter(fields...)
Fields:
cr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: CR time constant
RadiationDetectorDSP.CUSPChargeFilter
— Typestruct CUSPChargeFilter <: AbstractRadFIRFilter
CUSP filter.
For the definition the filter and a discussion of the filter properties, see
Constructors:
CUSPChargeFilter(; fields...)
Fields:
sigma::Union{Real, Unitful.AbstractQuantity{<:Real}}
: equivalent of shaping time (τₛ) Default: 450toplen::Union{Real, Unitful.AbstractQuantity{<:Real}}
: length of flat top (FT) Default: 10tau::Union{Real, Unitful.AbstractQuantity{<:Real}}
: decay constant of the exponential Default: 20length::Union{Real, Unitful.AbstractQuantity{<:Real}}
: total length of the filter (L) Default: 100beta::AbstractFloat
: scaling factor Default: 100.0
RadiationDetectorDSP.ConvolutionFilter
— Typestruct ConvolutionFilter{T<:RealQuantity} <: AbstractRadFIRFilter
A FIR filter defined by it's filter taps, applied via convolution with the input signal.
Constructors:
ConvolutionFilter(fields...)
Fields:
method::RadiationDetectorDSP.ConvolutionMethod
: Convolution methodcoeffs::AbstractVector{T} where T<:Union{Real, Unitful.AbstractQuantity{<:Real}}
: Filter tapsoffset::Int64
: Time axis offset
RadiationDetectorDSP.ConvolutionMethod
— Typeabstract type ConvolutionMethod
Indended as a type parameter to designate the behavior of a filter as linear or nonlinear.
Subtypes are DirectConvolution
and FFTConvolution
.
RadiationDetectorDSP.DNIMethod
— Typeabstract type DNIMethod
Abstract type for denoising and interpolation methods.
RadiationDetectorDSP.DifferentiatorFilter
— Typestruct DifferentiatorFilter <: AbstractRadIIRFilter
An integrator filter. It's inverse is IntegratorFilter
.
Constructors:
DifferentiatorFilter(fields...)
Fields:
gain::Union{Real, Unitful.AbstractQuantity{<:Real}}
: Filter gain
RadiationDetectorDSP.DirectConvolution
— TypeDirectConvolution() isa ConvolutionMethod
Compute filter convolutions directly, without FFT.
RadiationDetectorDSP.FFTConvolution
— TypeFFTConvolution() isa ConvolutionMethod
Compute filter convolutions via FFT.
RadiationDetectorDSP.FilteringType
— Typeabstract type FilteringType
Indended as a type parameter to designate the behavior of a filter as linear or nonlinear.
Subtypes are LinearFiltering
and NonlinearFiltering
.
RadiationDetectorDSP.FirstOrderIIR
— Typestruct FirstOrderIIR{T<:RealQuantity} <: AbstractRadIIRFilter
Constructors:
FirstOrderIIR(fields...)
Fields:
b_01::Tuple{T, T} where T<:Union{Real, Unitful.AbstractQuantity{<:Real}}
: Coefficients b0 to b1a_1::Tuple{T} where T<:Union{Real, Unitful.AbstractQuantity{<:Real}}
: Coefficient a1, a0 equals one implicitly
RadiationDetectorDSP.Gauss1DFilter
— Typestruct Gauss1DFilter <: AbstractRadFIRFilter
One dimensional gaussian filter defined as: f(x) = beta * exp(-0.5*(x/sigma)^2) / length
where x is in the interval [-alphasigma, alphasigma]
Constructors:
Gauss1DFilter(; fields...)
Fields:
sigma::Union{Real, Unitful.AbstractQuantity{<:Real}}
: standard deviation Default: 1.0length::Union{Real, Unitful.AbstractQuantity{<:Real}}
: total length of the filter Default: 100.0alpha::AbstractFloat
: the amount of standard deviations to cover in the gaussian window Default: 3.0beta::AbstractFloat
: scaling factor Default: 1.0
RadiationDetectorDSP.IntegratorCRFilter
— Typestruct IntegratorCRFilter <: AbstractRadIIRFilter
A modified CR-filter. The filter has an inverse.
Constructors:
IntegratorCRFilter(fields...)
Fields:
gain::Union{Real, Unitful.AbstractQuantity{<:Real}}
: Filter gaincr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: CR time constant
RadiationDetectorDSP.IntegratorFilter
— Typestruct IntegratorFilter <: AbstractRadIIRFilter
An integrator filter. It's inverse is DifferentiatorFilter
.
Constructors:
IntegratorFilter(fields...)
Fields:
gain::Union{Real, Unitful.AbstractQuantity{<:Real}}
: Filter gain
RadiationDetectorDSP.IntegratorModCRFilter
— Typestruct IntegratorModCRFilter <: AbstractRadIIRFilter
A modified CR-filter. The filter has an inverse.
Constructors:
IntegratorModCRFilter(fields...)
Fields:
gain::Union{Real, Unitful.AbstractQuantity{<:Real}}
: Filter gaincr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: CR time constant
RadiationDetectorDSP.Intersect
— Typestruct Intersect <: Function
Finds the intersects of a Y with a threshold
Constructors:
Intersect(; fields...)
Fields:
mintot::Union{Real, Unitful.AbstractQuantity{<:Real}}
: minimum time-over-threshold
RadiationDetectorDSP.InvCRFilter
— Typestruct InvCRFilter <: AbstractRadIIRFilter
Inverse of CRFilter
.
Constructors:
InvCRFilter(fields...)
Fields:
cr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: CR time constant
RadiationDetectorDSP.InvModCRFilter
— Typestruct InvModCRFilter <: AbstractRadIIRFilter
Inverse of ModCRFilter
.
Constructors:
InvModCRFilter(fields...)
Fields:
cr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: CR time constant
RadiationDetectorDSP.InvRCFilter
— Typestruct InvRCFilter <: AbstractRadIIRFilter
Inverse of RCFilter
.
Constructors:
InvRCFilter(fields...)
Fields:
rc::Union{Real, Unitful.AbstractQuantity{<:Real}}
: RC time constant
RadiationDetectorDSP.InvSecondOrderCRFilter
— Typestruct InvSecondOrderCRFilter <: AbstractRadIIRFilter
Inverse of SecondOrderCRFilter
. Apply a double pole-zero cancellation using the provided time constants to the waveform.
Constructors:
InvSecondOrderCRFilter(fields...)
Fields:
cr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: time constant of the first exponential to be deconvolvedcr2::Union{Real, Unitful.AbstractQuantity{<:Real}}
: time constant of the second exponential to be deconvolvedf::Real
: the fraction faktor which the second exponential contributes
RadiationDetectorDSP.LinearFiltering
— Typeabstract type LinearFiltering <: FilteringType
When used as a type parameter value, marks linear behavior of a filter.
RadiationDetectorDSP.ModCRFilter
— Typestruct ModCRFilter <: AbstractRadIIRFilter
A first-order CR highpass filter, modified for full-amplitude step-signal response.
The resonse of the standard digital CRFilter
will not recover the full amplitude of a digital step stignal since a step from one sample to the still has a finite rise time. This version of a CR filter compensates for this loss in amplitude, so it effectively treats a step as having
The inverse filter is InvModCRFilter
, this is typically stable even in the presence of additional noise (see CRFilter
).
Constructors:
ModCRFilter(fields...)
Fields:
cr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: CR time constant
RadiationDetectorDSP.NonlinearFiltering
— Typeabstract type NonlinearFiltering <: FilteringType
When used as a type parameter value, marks non linear behavior of a filter.
RadiationDetectorDSP.PolynomialDNI
— Typestruct PolynomialDNI <: DNIMethod
Polynomial denoising and interpolation method.
Operates in a similar way as a Savitzky-Golay filter, but interpolates as well.
Constructors:
PolynomialDNI(; fields...)
Fields:
degree::Int64
: polynomial degree Default: (2, "length")length::Union{Real, Unitful.AbstractQuantity{<:Real}}
RadiationDetectorDSP.RCFilter
— Typestruct RCFilter <: AbstractRadII>RFilter
A first-order RC lowpass filter.
The inverse filter is InvCRFilter
, but note that this is unstable in the presence of additional noise. As an RC filter attenuates high-frequency noise, its inverse amplifies such noise and will typically not be useful to deconvolve signals in practical applications.
Constructors:
RCFilter(fields...)
Fields:
rc::Union{Real, Unitful.AbstractQuantity{<:Real}}
: RC time constant
RadiationDetectorDSP.SamplingInfo
— Typestruct SamplingInfo{T<:RealQuantity,A<:AbstractVector{<:RealQuantity}}
Holds sampling information.
The numerical type of an individual sample is T
, the (time) axis is given by the axis
field.
Constructors:
SamplingInfo{T,A}(axis)
SamplingInfo{T}(axis)
Fields:
axis::AbstractVector{<:Union{Real, Unitful.AbstractQuantity{<:Real}}}
RadiationDetectorDSP.SavitzkyGolayFilter
— Typestruct SavitzkyGolayFilter <: AbstractRadFIRFilter
Constructors:
SavitzkyGolayFilter(; fields...)
Fields:
length::Union{Real, Unitful.AbstractQuantity{<:Real}}
: filter lengthdegree::Int64
: Polynomial defgreederivative::Int64
: n-th derivative (0 for no derivative)
RadiationDetectorDSP.SecondOrderCRFilter
— Typestruct SecondOrderCRFilter <: AbstractRadIIRFilter
A scond order CR highpass filter. The filter has an inverse InvSecondOrderCRFilter
.
Constructors:
SecondOrderCRFilter(fields...)
Fields:
cr::Union{Real, Unitful.AbstractQuantity{<:Real}}
: time constant of the first exponential to be deconvolvedcr2::Union{Real, Unitful.AbstractQuantity{<:Real}}
: time constant of the second exponential to be deconvolvedf::Real
: the fraction faktor which the second exponential contributes
RadiationDetectorDSP.SignalEstimator
— Typestruct SignalEstimator <: Function
Estimates a signal at a given position x
.
Usage:
(f::SamplesOrWaveform)(input::RDWaveform, x::RealQuantity)
Constructors:
SignalEstimator(; fields...)
Fields:
dni::DNIMethod
: denoising and interpolation method
RadiationDetectorDSP.SimpleCSAFilter
— Typestruct SimpleCSAFilter <: AbstractRadIIRFilter
Simulates the current-signal response of a charge-sensitive preamplifier with resistive reset, the output is a charge signal.
It is equivalent to the composition
CRFilter(cr = tau_decay) ∘
Integrator(gain = gain) ∘
RCFilter(rc = tau_rise)
and maps to a single BiquadFilter
.
This filter has an inverse, but the inverse is very unstable in the presence of additional noise if tau_rise
is not zero (since the inverse of an RC-filter is unstable under noise). Even if tau_rise
is zero the inverse will still amplify noise (since it differentiates), so it should be used very carefully when deconvolving signals in practical applications.
Constructors:
SimpleCSAFilter(fields...)
Fields:
tau_rise::Union{Real, Unitful.AbstractQuantity{<:Real}}
: Rise time constanttau_decay::Union{Real, Unitful.AbstractQuantity{<:Real}}
: Decay time constantgain::Union{Real, Unitful.AbstractQuantity{<:Real}}
: Gain
RadiationDetectorDSP.TrapezoidalChargeFilter
— Typestruct TrapezoidalChargeFilter <: AbstractRadNonlinearFilter
Filter that responds to a step signal with a trapezoidal pulse.
The filter is equivalent to two moving averages separated by a gap.
Constructors:
TrapezoidalChargeFilter(; fields...)
Fields:
avgtime::Union{Real, Unitful.AbstractQuantity{<:Real}}
: pre-rise averaging timegaptime::Union{Real, Unitful.AbstractQuantity{<:Real}}
: gap timeavgtime2::Union{Real, Unitful.AbstractQuantity{<:Real}}
: post-rise averaging time
A sharp step on the input will result in a trapezoid with rise time and fall time avgtime
and a flat top of length gaptime
.
RadiationDetectorDSP.TruncateFilter
— Typestruct TruncateFilter <: AbstractRadSigFilter{LinearFiltering}
Filter that truncates the input signal.
Constructors:
TruncateFilter(; fields...)
Fields:
interval::IntervalSets.ClosedInterval{<:Union{Real, Unitful.AbstractQuantity{<:Real}}}
: interval to keep
RadiationDetectorDSP.ZACChargeFilter
— Typestruct ZACChargeFilter <: AbstractRadFIRFilter
Zero area cusp (ZAC) filter.
For the definition the filter and a discussion of the filter properties, see "Improvement of the energy resolution via an optimized digital signal processing in GERDA Phase I", Eur. Phys. J. C 75, 255 (2015).
Constructors:
ZACChargeFilter(; fields...)
Fields:
sigma::Union{Real, Unitful.AbstractQuantity{<:Real}}
: equivalent of shaping time (τₛ) Default: 450toplen::Union{Real, Unitful.AbstractQuantity{<:Real}}
: length of flat top (FT) Default: 10tau::Union{Real, Unitful.AbstractQuantity{<:Real}}
: decay constant of the exponential Default: 20length::Union{Real, Unitful.AbstractQuantity{<:Real}}
: total length of the filter (L) Default: 100beta::AbstractFloat
: scaling factor Default: 100.0
RadiationDetectorDSP.SamplesOrWaveform
— Typeconst RadiationDetectorDSP.SamplesOrWaveform{T<:RealQuantity} = Union{AbstractSamples{T},RDWaveform{<:Any,T}}
A vector of signal samples or a waveform.
RadiationDetectorDSP.adapt_memlayout
— FunctionRadiationDetectorDSP::adapt_memlayout(
fi::AbstractRadSigFilterInstance,
backend::KernelAbstractions.Backend,
A::AbstractArray{<:Number}
)
Adapts the memory layout of A
in a suitable fashion for fi
on computing device backend
.
Returns a row-major version of A
on all backends by default, filter instance types may specialize this behavior.
RadiationDetectorDSP.add_rect_pulse!
— Functionadd_rect_pulse!(samples::AbstractSamples, start::Integer, pulselen::Integer, amplitude::Real = 1.0)
Add a rectangular pulse to samples
.
RadiationDetectorDSP.bc_rdfilt
— Functionbc_rdfilt(flt::AbstractRadSigFilter, input)
bc_rdfilt(fi::AbstractRadSigFilterInstance, input)
Broadcast filter instance fi
over signals input
, return the filtered signals.
RadiationDetectorDSP.bc_rdfilt!
— Functionbc_rdfilt!(outputs, fi::AbstractRadSigFilterInstance, inputs)
Broadcast filter flt
or filter instance fi
over signals inputs
, storing the results in outputs
.
inputs
and outputs
must be of type AbstractVector{<:AbstractSamples}
.
Returns outputs
.
RadiationDetectorDSP.charge_trapflt!
— Methodcharge_trapflt!(samples::AbstractVector{<:RealOrSIMD{<:AbstractFloat}}, navg::Integer, ngap::Integer)
Apply a trapezoidal FIR filter to a charge signal in samples
.
RadiationDetectorDSP.cr_filter
— Methodcr_filter(CR::Real)
Return a DSP.jl-compatible CR-filter.
RadiationDetectorDSP.crmod_filter
— Methodcrmod_filter(CR::Real)
Return a DSP.jl-compatible modified CR-filter.
RadiationDetectorDSP.cusp_charge_filter_coeffs
— Methodcusp_charge_filter_coeffs(N::Int, sigma::U, FT::Int, tau::V, beta::W
) where {U, V, W <: AbstractFloat}
return a vector representing the cusp filter applicaible on a charge signal, where N
is the total length of the filter, FT
the length of the flat top, sigma
the filter shaping time,tau
the decay constant and a
the scaling factor.
RadiationDetectorDSP.dfilt!
— Functionrdfilt!(output, fi::AbstractRadSigFilterInstance, input)
Apply filter flt
or filter instance fi
to signal input
and store the filtered signal in output
. Return output
.
RadiationDetectorDSP.differentiator_filter
— Methoddifferentiator_filter(gain::Real)
Return a DSP.jl-compatible differentiator filter.
RadiationDetectorDSP.elsmplinfo
— Functionsmplinfo(smpls::AbstractSamples)::SamplingInfo
smplinfo(wf::RDWaveform{T,U})::RDWaveform
Get sampling information an array of vectors of samples, resp. an array of waveform. All elements must have equal samling information.
RadiationDetectorDSP.flt_input_length
— FunctionRadiationDetectorDSP.flt_input_length(fi::AbstractRadSigFilterInstance)::Integer
Get the output signal length of a filter instance fi
.
Must be implemented for all subtypes of AbstractRadSigFilterInstance
.
RadiationDetectorDSP.flt_input_smpltype
— FunctionRadiationDetectorDSP.flt_input_smpltype(fi::AbstractRadSigFilterInstance)
Get the input sample type of a filter instance fi
.
Must be implemented for all subtypes of AbstractRadSigFilterInstance
.
RadiationDetectorDSP.flt_output_length
— FunctionRadiationDetectorDSP.flt_output_length(fi::SomeFilterInstance)::Integer
Get the output signal length of filter instance fi
.
Must be implemented for all subtypes of AbstractRadSigFilterInstance
.
RadiationDetectorDSP.flt_output_smpltype
— FunctionRadiationDetectorDSP.flt_output_smpltype(fi::AbstractRadSigFilterInstance)
Get the output sample type for
- a filter
flt
given an input sample typeinput_smpltype
- a filter instance
fi
Must be implemented for all subtypes of AbstractRadSigFilter
.
RadiationDetectorDSP.flt_output_time_axis
— FunctionRadiationDetectorDSP.flt_output_time_axis(fi::SomeFilterInstance, time::AbstractVector{<:RealQuantity})::AbstractVector{<:RealQuantity}
Get the output time axis of a filter instance fi
, given an input time axis time
.
Must be implemented for subtypes of AbstractRadSigFilter
and AbstractRadSigFilterInstance
only if the filter's output time axis can be computed directly from the input time axis.
RadiationDetectorDSP.fltinstance
— Functionfltinstance(flt::AbstractRadSigFilter, si::SamplingInfo)::AbstractRadSigFilterInstance
Create a filter instance of the filter flt
, specialized for the given input, resp. input characteristics.
RadiationDetectorDSP.gaussian_coeffs
— Methodgaussian_coeffs(N::Int, sigma::V, alpha::U, beta::U) where {V, U}
compute a gaussian kernel, where N
is the total length of the kernel, alpha
the amount of standard deviations to cover, sigma
the standard deviation and beta
the total scaling factor.
RadiationDetectorDSP.gen_rect_pulse
— Functiongen_rect_pulse(tracelen::Integer, start::Integer, pulselen::Integer, amplitude::Real = 1.0)
Generate a rectangular pulse.
RadiationDetectorDSP.integrator_cr_filter
— Methodintegrator_cr_filter(gain::Real, CR::Real)
Return a DSP.jl-compatible integrator plus CR filter.
RadiationDetectorDSP.integrator_crmod_filter
— Methodintegrator_crmod_filter(gain::Real, CR::Real)
Return a DSP.jl-compatible integrator plus modified CR filter.
RadiationDetectorDSP.integrator_filter
— Methodintegrator_filter(gain::Real)
Return a DSP.jl-compatible integrator filter.
RadiationDetectorDSP.inv_cr_filter
— Methodinv_cr_filter(CR::Real)
Return a DSP.jl-compatible inverse CR-filter.
RadiationDetectorDSP.inv_crmod_filter
— Methodinv_crmod_filter(CR::Real)
Return a DSP.jl-compatible inverse modified CR-filter.
RadiationDetectorDSP.inv_rc_filter
— Methodinv_rc_filter(RC::Real)
Return a DSP.jl-compatible RC-filter.
RadiationDetectorDSP.multiply_waveform
— Functionmultiply_waveform(signal::AbstractSamples, a::RealQuantity)
multiply_waveform(wf::RDWaveform, a::RealQuantity)
Multiplies each sample of a waveform by a
.
RadiationDetectorDSP.rc_filter
— Methodrc_filter(RC::Real)
Return a DSP.jl-compatible RC-filter.
RadiationDetectorDSP.rdfilt
— Functionrdfilt(fi::AbstractRadSigFilterInstance, input)
Apply filter instance fi
to signal input
, return the filtered signal.
Returns output
.
RadiationDetectorDSP.reverse_waveform
— Functionreverse_waveform(signal::AbstractSamples)
reverse_waveform(wf::RDWaveform)
Reverses the order of samples in a waveform.
RadiationDetectorDSP.shift_waveform
— Functionshift_waveform(signal::AbstractSamples, a::RealQuantity)
shift_waveform(wf::RDWaveform, a::RealQuantity)
Shifts each sample of a waveform up by a
.
RadiationDetectorDSP.simple_csa_response_filter
— Functionsimplecsaresponsefilter(τrise::Real, τdecay::Real, gain::Real = one(τrise))
Return a DSP.jl-compatible filter that models the response of a typical charge-sensitive amplifier (CSA).
RadiationDetectorDSP.smplinfo
— Functionsmplinfo(smpls::AbstractSamples)::SamplingInfo
smplinfo(wf::RDWaveform{T,U})::RDWaveform
Get sampling information from a vector of samples, resp. a waveform.
RadiationDetectorDSP.zac_charge_filter_coeffs
— Methodzac_charge_filter_coeffs(
N::Int, sigma::V, FT::Int, tau::T, beta::U
) where {V, T, U <: AbstractFloat}
return a vector representing the zac filter applicaible on a charge signal, where N
is the total length of the filter, FT
the length of the flat top, sigma
the filter shaping time, tau
the decay constant and beta an additional scaling factor. (see Eur. Phys. J. C (2015) 75:255).