API
Measurement data processing
AntibodyMethodsDoseResponse.dose_response_check
— Functiondose_response_check(fitting_data::FittingData)
Test if the FittingData
object fitting_data
satisfies properties for dose-response data:
- only real numbers
- numbers must be positive
- no
NaN
orInf
in the data set
AntibodyMethodsDoseResponse.normalize_data
— Functionnormalize_data!(fitting_data::FittingData; offset::Real = 0, reference::Union{Nothing,Real} = nothing)
Return the FittingData
object fitting_data
with normalized response data and errors.
reference = nothing
: The reference signal to normalize the responses to. Ifreference = nothing
, the maximal response is used a reference point.offset = 0
: Signal offset to be subtracted (applies to both the responses and the reference point).
The normalization of the responses
and errors
are:
new_responses = (responses - offset) / (reference - offset)
errors = errors / (reference - offset)
AntibodyMethodsDoseResponse.normalize_data!
— Functionnormalize_data!(fitting_data::FittingData; offset::Real = 0, reference::Union{Nothing,Real} = nothing)
Normalize the FittingData
object fitting_data
by mutation.
reference = nothing
: The reference signal to normalize the responses to. Ifreference = nothing
, the maximal response is used a reference point.offset = 0
: Signal offset to be subtracted (applies to both the responses and the reference point).
The normalization of the responses
and errors
are:
new_responses = (responses - offset) / (reference - offset)
errors = errors / (reference - offset)
Depletion correction of measurement data
AntibodyMethodsDoseResponse.scale_bound
— Functionscale_bound(fitting_data::FittingData)
Get the upper bound for the scale factor β
s.t. a_i - β*r_i ≥ 0
, where r_i
are the responses and a_i
are the (initial) antibody concentrations.
AntibodyMethodsDoseResponse.simple_depletion_correction
— Functionsimple_depletion_correction(fitting_data::FittingData,scale::Real)
Return the depletion-corrected FittingData
object.
The concentrations are corrected to concentration - scale * response
.
simple_depletion_correction(fitting_data::FittingData)
Return depletion-corrected FittingData
object using the largest possible scale factor.
The concentrations are corrected to concentration - scale * response
.
AntibodyMethodsDoseResponse.simple_depletion_correction!
— Functionsimple_depletion_correction!(fitting_data::FittingData,scale::Real)
Depletion-correct the FittingData
object fitting_data
by mutation.
The concentrations are corrected to concentration - scale * response
.
simple_depletion_correction!(fitting_data::FittingData)
Depletion-correct the FittingData
object fitting_data
by mutation, using the largest possible scale factor.
The concentrations are corrected to concentration - scale * response
.
Models
AntibodyMethodsDoseResponse.accumulation_model
— Functionaccumulation_model(grid::OneDimGrid; offset = nothing)
Create a multi-epitope accumulation model. Returns (model,λ,centers,volumes)
where
model
is aModelFunctions
object.λ
is an initial parameter array (the weights of the grid and the offset ifoffset != nothing
). Ifoffset != nothing
, the last element is the offset parameterλ[end] = offset
.centers
andvolumes
are the remaining properties of the grid, seeexport_all
Model function
The following model function and partial derivatives are used:
\[\text{model}(a,\lambda) = \lambda_e + \sum_i \lambda_i \left(1-e^{-\frac{a}{c_i}}\right) \approx \lambda_e + \sum_i \int_{l_i}^{u_i} \frac{\lambda_i}{u_i-l_i}\left(1-e^{-\frac{a}{k}}\right) \ dk \]
\[\partial_{\lambda_j} \text{model}(a,\lambda) = 1-e^{-\frac{a}{c_i}} \quad, \qquad \partial_{\lambda_e} \text{model}(a,\lambda) = 1\]
where $a$ is the antibody concentration, $c_i$ are the centers of the grid
intervals $[u_i,l_i]$ and $\lambda_e$ is the offset (if offset != nothing
).
AntibodyMethodsDoseResponse.langmuir_model
— Functionlangmuir_model(grid::OneDimGrid; offset = nothing)
Create a multi-epitope Langmuir model. Returns (model,λ,centers,volumes)
where
model
is aModelFunctions
object.λ
is an initial parameter array (the weights of the grid and the offset ifoffset != nothing
). Ifoffset != nothing
, the last element is the offset parameterλ[end] = offset
.centers
andvolumes
are the remaining properties of the grid, seeexport_all
.
Model function
The following model function and partial derivatives are used:
\[\text{model}(a,\lambda) = \lambda_e + \sum_i \frac{\lambda_i\cdot a}{ (u_i - l_i)} \ln\left(\frac{a+u_i}{a+l_i}\right) = \lambda_e + \sum_i \int_{l_i}^{u_i} \frac{\frac{\lambda_i}{u_i - l_i}}{1+\frac{k}{a}} \ dk \]
\[\partial_{\lambda_j} \text{model}(a,\lambda) = \frac{a}{u_j-l_j} \ln\left(\frac{a+u_j}{a+l_j}\right)\quad, \qquad \partial_{\lambda_e} \text{model}(a,\lambda) = 1\]
where $a$ is the antibody concentration, $[l_i,u_i]$ are the intervals of the grid
and $\lambda_e$ is the offset (if offset != nothing
).
AntibodyMethodsDoseResponse.accumulation_inv_const_model
— Functionaccumulation_inv_const_model(grid::OneDimGrid; offset = nothing)
Create a multi-epitope accumulation model with 1/K_τ = k_a * τ
as constant domain. Returns (model,λ,centers,volumes)
where
model
is aModelFunctions
object.λ
is an initial parameter array (the weights of the grid and the offset ifoffset != nothing
). Ifoffset != nothing
, the last element is the offset parameterλ[end] = offset
.centers
andvolumes
are the remaining properties of the grid, seeexport_all
Model function
The following model function and partial derivatives are used:
\[\text{model}(a,\lambda) = \lambda_e + \sum_i \lambda_i \left(1+\frac{1}{a\cdot(u_i-l_i)}\left(e^{-a u_i}-e^{-a l_i} \right) \right) = \lambda_e + \sum_i \int_{l_i}^{u_i} \frac{\lambda_i}{u_i-l_i}\left(1-e^{-a k}\right) \ dk \]
\[\partial_{\lambda_j} \text{model}(a,\lambda) = 1+\frac{1}{a\cdot(u_i-l_i)}\left(e^{-a u_i}-e^{-a l_i} \right) \quad, \qquad \partial_{\lambda_e} \text{model}(a,\lambda) = 1\]
where $a$ is the antibody concentration, $c_i$ are the centers of the grid
intervals $[u_i,l_i]$ and $\lambda_e$ is the offset (if offset != nothing
).
AntibodyMethodsDoseResponse.langmuir_inv_const_model
— Functionlangmuir_inv_const_model(grid::OneDimGrid; offset = nothing)
Create a multi-epitope Langmuir model with 1/K_d = k_a / k_d
as constant domain. Returns (model,λ,centers,volumes)
where
model
is aModelFunctions
object.λ
is an initial parameter array (the weights of the grid and the offset ifoffset != nothing
). Ifoffset != nothing
, the last element is the offset parameterλ[end] = offset
.centers
andvolumes
are the remaining properties of the grid, seeexport_all
Model function
The following model function and partial derivatives are used:
\[\text{model}(a,\lambda) = \lambda_e + \sum_i \lambda_i \left(1+ \frac{1}{a\cdot (u_i-l_i)} \ln\left(\frac{a l_i +1}{a u_i +1}\right)\right) = \lambda_e + \sum_i \int_{l_i}^{u_i} \frac{\frac{\lambda_i}{u_i - l_i}}{1+\frac{1}{a\cdot k}} \ dk \]
\[\partial_{\lambda_j} \text{model}(a,\lambda) = 1+ \frac{1}{a\cdot (u_i-l_i)} \ln\left(\frac{a l_i +1}{a u_i +1}\right)\quad, \qquad \partial_{\lambda_e} \text{model}(a,\lambda) = 1\]
where $a$ is the antibody concentration, $[l_i,u_i]$ are the intervals in the grid
and $\lambda_e$ is the offset (if offset != nothing
).
Result type
AntibodyMethodsDoseResponse.DoseResponseResult
— Typestruct DoseResponseResult
Data type to store dose-response result data (e.g. from a dose-response-curve fitting or a simulation).
Fields
concentrations
: Antibody concentrations.responses
: Response values.
Default constructor
DoseResponseResult(concentrations,responses)
Model constructor
DoseResponseResult(grid::OneDimGrid,concentrations;
offset::Real = 0,
model::Function = accumulation_model
)
Calculate a dose-response curve from a K_τ grid
and a model
for given concentrations
. The offset
value is a global additive shift for all response values.
The available model functions are accumulation_model
, accumulation_inv_const_model
, langmuir_model
and langmuir_inv_const_model
.
Adaptive fitting
AntibodyMethodsDoseResponse.AdaptiveOptions
— Typemutable struct AdaptiveOptions
Data type to define adaptive_dose_response_fit
options.
Constructor
AdaptiveOptions(keywords...)
The following keywords (with default values) are available:
name::AbstractString = "Adaptive optimization"
: The name that is used whenshow_progress==true
.show_progress::Bool = true
: Show progress in standard output.iterations::Integer = 1
: Number of refinement iterations.model::Function = accumulation_model
: The model-function that is used for the data-fit. The available model functions areaccumulation_model
,accumulation_inv_const_model
,langmuir_model
andlangmuir_inv_const_model
.offset = nothing
: Offset parameter for the model function. Ifnothing
, no offset is used.objective::Symbol = :lsq
. The objective function for the data-fit. Available are:lsq
,:posterior
and:log_posterior
.prior_generator::Function = default_prior_generator
: The function that generates the prior. The function must have the signature(grid_centers,grid_volumes,offset)
and must return a functionλ-> prior(λ)
orλ-> log_prior(λ)
in case of a:log_posterior
objective. Thedefault_prior_generator
generates a uniform priorλ-> 0
for the log-posterior objective.distribution_derivatives = nothing
: Array of partial derivatives of the logarithmic distributions for the log-posterior objective. Seelog_posterior_gradient
.prior_gradient_generator = default_prior_gradient_generator
: The function that generates the log-prior gradient (seelog_posterior_gradient
). The function must have the signature(grid_centers,grid_volumes,offset)
and must return a functionλ-> ∇log_prior(λ)
. Thedefault_prior_gradient_generator
returnsnothing
which internally corresponds to the uniform prior for the log-posterior objective.block_variation::Function =
log_area_scaled_variation
andselection::Function = maximum
are the refinement options ofrefine!
.
AntibodyMethodsDoseResponse.area_scaled_variation
— Functionarea_scaled_variation(center, volume, weight,
neighbor_centers, neighbor_volumes, neighbor_weights)
Block variation function for refine!
. Variation value based on the difference of the weights, scaled with the area (volume) of the corresponding block.
mean(@. abs(weight * volume - neighbor_weights * neighbor_volumes))
AntibodyMethodsDoseResponse.log_area_scaled_variation
— Functionlog_area_scaled_variation(center, volume, weight,
neighbor_centers, neighbor_volumes, neighbor_weights)
Block variation function for refine!
. Variation value based on the difference of the weight, scaled with the visible area (in a logarithmic plot) of the corresponding block.
log_volume = (log10(center + volume / 2) - log10(center - volume / 2))
neighbor_log_volumes = @. (log10(neighbor_centers + neighbor_volumes / 2) - log10(neighbor_centers - neighbor_log_volumes / 2))
mean(@. abs(weight * log_volume - neighbor_weights * neighbor_log_volumes))
AntibodyMethodsDoseResponse.adaptive_dose_response_fit
— Functionadaptive_dose_response_fit(initial_grid::OneDimGrid,
data::FittingData,
minimizer::Function;
options::AdaptiveOptions=AdaptiveOptions()
)
Fit dose-response data
(with adaptive grid refinements depending on the options
) and return an AdaptiveResult
object. The initial gird is not mutated. For the AdaptiveResult
object a copy of the gird is created.
Minimizer function
- The sign of posterior and log-posterior objectives is flipped for consistency reasons. A minimizer needs to be used for all objectives (
:lsq
,:posterior
,:log_posterior
). minimizer
: The function that minimizes the objective function. It needs to be specified by the user and must have the signature(objective_function, objective_gradient!,parameters)
.- The
objective_function
always has the signature(parameters)
. - The
objective_gradient!
can benothing
. Otherwise it must have the signature(gradient_vector, parameter)
. It must mutate thegradient_vector
and return the mutatedgradient_vector
. parameters
is the initial parameter array to start minimization from.
Gradients for minimization
Whether objective_gradient!
is nothing
or a proper function depends on the specified options
(see AdaptiveOptions
).
- The
:lsq
objective always provides analytical gradients. - The
:posterior
objective never provides analytical gradients. - The
:log_posterior
objective only provides analytical gradients ifdistribution_derivatives != nothing
.
AntibodyMethodsDoseResponse.AdaptiveResult
— Typemutable struct AdaptiveResult
Data type used by adaptive_dose_response_fit
to summarize the results.
The struct has the following fields:
result
: TheDoseResponseResult
object corresponding to the fit result.grid
: The grid (with imported weights) corresponding to the fit result.optimizer
: The raw result parameter.objective_value
: The objective-function value ofoptimizer
.time
: The elapsed time for the model fit (in seconds).
Convenience methods
AntibodyMethodsDoseResponse.LogRange
— FunctionLogRange(start::Real,stop::Real,n::Integer, base::Real = 10.0)
Return a vector with n
values logarithmically distributed between start
and stop
. The logarithm base
can be changed with the last, optional argument.
AntibodyMethodsDoseResponse.peak_detection
— Functionpeak_detection(grid::OneDimGrid, relative_threshold::Real = 0.1;
volume_normalization::Symbol = :log,
fill::Bool = true
)
Return (peak_group_indices,peak_group_domains)
w.r.t. the relativ_ threshold
.
peak_group_indices
contains the index-vectors of the peaks. The index order is the order ofexport_weights
andexport_all
.peak_group_domains
contains the intervals covered by the respective peaks.- The cutoff threshold is determined by
relative_threshold * largest_weight
. - If
fill == true
, the gaps between the peaks are added topeak_group_indices
andpeak_group_domains
. volume_normalization
normalizes the weights of the grid (without mutation) before the peaks are determined.:none
uses the raw weights.:linear
divides the weights by the block volume and:log
divides the weight by the block volume in a logarithmic scale.
Uncertainty estimation
AntibodyMethodsDoseResponse.EpitopeUncertainty
— Typestruct EpitopeUncertainty
Data type to store uncertainty estimates for the weights
of a K_τ grid
.
Fields
levels
: List of uncertainty levels. They can be sample quantiles or fractions of the best objective value, depending on the constructor.lower
: Matrix of estimated lower bounds for the weights at the corresponding uncertainty levels (dimension order : [level, grid parameter index]).upper
: Matrix of estimated upper bounds for the wights at the corresponding uncertainty levels (dimension order : [level, gird parameter index]).lower_offset
: Vector of lower bounds for the offset parameter corresponding to the uncertainty levels.upper_offset
: Vector of upper bounds for the offset parameter corresponding to the uncertainty levels.
Default constructor
EpitopeUncertainty(levels, lower, upper, lower_offset= nothing, upper_offset = nothing)
Construction from bin-wise shifting
EpitopeUncertainty(data::FittingData,grid::OneDimGrid, bins = collect(1:length(grid));
keywords...)
Estimate uncertainty by shifting all grid weights uniformly, one bin at a time, while keeping the other bins fixed. Admissible weights (for a given level) are determined by calculating the objective-function value (objective function automatically generated) for the shifted weights.
The bins can be defined as vector of indices, e.g. [[1,2,3],[4,5,6]]
or [1,3,5]
which is converted to [[1],[3],[5]]
. To obtain bin indices from the grid-domain, use select_indices
.
The following keywords are available:
levels = collect(0.1:0.1:1)
: The uncertainty levels as fractions of the best objective value.steps::Integer = 10^4
: Number of intermediate shifts to be tested. The maximal shift range is determined automatically.bisections::Integer = 10^2
: Number of interval-bisection steps to determine the maximal shift range.volume_normalization = :none
: Changes the scaling of the weight shifts for each individual weight in a bin. Use:none
to apply the same shift for each weight. Use:linear
to scale the shifts with the interval volumes corresponding to the weights. And use:log
to scale the sifts with the visual interval volumes (as they appear in a logarithmic plot) corresponding to the weights.options::AdaptiveOptions = AdaptiveOptions()
: The objective function is automatically generated, using the same construction asadaptive_dose_response_fit
. SeeAdaptiveOptions
andadaptive_dose_response_fit
for the details. To use the offset estimation from a fit result, add the estimated value to theAdaptiveOptions
with theoffset
keyword:options = AdaptiveOptions(other_options..., offset = estimated_offset)
.
Construction from samples
EpitopeUncertainty(samples; keywords... )
Estimate the uncertainty as credibility intervals (symmetric interval around the median) from samples
(drawn from a posterior distribution). The samples can either be passed as array of samples (parameters) or as matrix (order: [parameter index, sample index]).
The following keywords are available:
levels = collect(0.1:0.1:1)
: The uncertainty levels as quantiles.offset::Bool = false
: Iftrue
the last parameter element is treated as offset parameter, otherwise, all parameters are treated as grid weights.
AntibodyMethodsDoseResponse.DoseResponseUncertainty
— Typestruct DoseResponseUncertainty
Data type to store uncertainty estimates for the response values of a DoseResponseResult
object.
Fields
levels
: List of uncertainty levels, corresponding to the uncertainty levels of aEpitopeUncertainty
object (if constructed from anEpitopeUncertainty
object).concentration
: Vector of concentrations (no uncertainty).lower
: Matrix of estimated lower bounds for the responses at the corresponding uncertainty level (dimension order : [level, concentration/response index]).upper
: Matrix of estimated upper bounds for the responses at the corresponding uncertainty level (dimension order : [level, concentration/response index]).
Default constructor
DoseResponseUncertainty(levels,concentrations,lower,upper)
Construction from an EpitopeUncertainty object
DoseResponseUncertainty(grid::OneDimGrid,
eu::EpitopeUncertainty,
concentrations::AbstractVector;
keywords...
)
Estimate the dose-response uncertainty from an EpitopeUncertainty
object eu
for the provided concentrations
. The grid
should be the grid that was used to create the EpitopeUncertainty
object eu
.
The following keywords are available:
bins = [collect(1:length(grid))]
: The response bounds are calculated as point-wise minima/maxima of responses created from the grid weights, where one bin at a time is replaced with theEpitopeUncertainty
lower and upper bound, while keeping the other weights fixed. For the minima/maxima all response values, iterating over all bins, are considered. Ideally, the bins should correspond to the bins that were used to construct theEpitopeUncertainty
objecteu
.model::Function = accumulation_model
: The model that is used to calculate the response values. The available model functions areaccumulation_model
,accumulation_inv_const_model
,langmuir_model
andlangmuir_inv_const_model
.
There is no offset
keyword, as the offsets are determined by the EpitopeUncertainty
object.