Docstrings
Here we give some of the main docstrings.
LikelihoodProblem
ProfileLikelihood.AbstractLikelihoodProblem
— Typeabstract type AbstractLikelihoodProblem{N, L}
Abstract type of a likelihood problem, where N
is the number of parameters and L
is the type of the likelihood function.
ProfileLikelihood.LikelihoodProblem
— TypeLikelihoodProblem{N,P,D,L,Θ,S} <: AbstractLikelihoodProblem
Struct representing a likelihood problem.
Fields
problem::P
The associated OptimizationProblem
.
data::D
The argument p
used in the log-likelihood function.
log_likelihood_function::L
The log-likelihood function, taking the form ℓ(θ, p)
.
θ₀::Θ
Initial estimates for the MLE θ
.
syms::S
Variable names for the parameters. This is wrapped into a SymbolCache
from SymbolicIndexingInterface.jl.
The extra parameter N
is the number of parameters.
Constructors
Standard
LikelihoodProblem(loglik::Function, θ₀;
syms=eachindex(θ₀), data=SciMLBase.NullParameters(),
f_kwargs=nothing, prob_kwargs=nothing)
Constructor for the LikelihoodProblem
.
Arguments
loglik::Function
: The log-likelihood function, taking the formℓ(θ, p)
.θ₀
: The estimates estimates for the MLEs.
Keyword Arguments
syms=eachindex(θ₀)
: Names for each parameter.data=SciMLBase.NullParameters()
: The parameterp
in the log-likelihood function.f_kwargs=nothing
: Keyword arguments, passed as aNamedTuple
, for theOptimizationFunction
.prob_kwargs=nothing
: Keyword arguments, passed as aNamedTuple
, for theOptimizationProblem
.
Outputs
Returns the LikelihoodProblem
problem object.
With arguments for a differential equation problem
LikelihoodProblem(loglik::Function, θ₀,
ode_function, u₀, tspan;
syms=eachindex(θ₀), data=SciMLBase.NullParameters(),
ode_parameters=SciMLBase.NullParameters(), ode_alg,
ode_kwargs=nothing, f_kwargs=nothing, prob_kwargs=nothing)
Constructor for the LikelihoodProblem
for a differential equation problem.
Arguments
loglik::Function
: The log-likelihood function, taking the formℓ(θ, p, integrator)
.θ₀
: The estimates estimates for the MLEs.ode_function
: The functionf(du, u, p, t)
orf(u, p, t)
for the differential equation.u₀
: The initial condition for the differential equation.tspan
: The time-span to solve the differential equation over.
Keyword Arguments
syms=eachindex(θ₀)
: Names for each parameter.data=SciMLBase.NullParameters()
: The parameterp
in the log-likelihood function.ode_parameters=SciMLBase.NullParameters()
: The parameterp
inode_function
.ode_alg
: The algorithm used for solving the differential equatios.ode_kwargs=nothing
: Extra keyword arguments, passed as aNamedTuple
, to pass into the integrator; seeconstruct_integrator
.f_kwargs=nothing
: Keyword arguments, passed as aNamedTuple
, for theOptimizationFunction
.prob_kwargs=nothing
: Keyword arguments, passed as aNamedTuple
, for theOptimizationProblem
.
Outputs
Returns the LikelihoodProblem
problem object.
With an integrator
LikelihoodProblem(loglik::Function, θ₀, integrator;
syms=eachindex(θ₀), data=SciMLBase.NullParameters(),
f_kwargs=nothing, prob_kwargs=nothing)
Constructor for the LikelihoodProblem
for a differential equation problem with associated integrator
.
Arguments
loglik::Function
: The log-likelihood function, taking the formℓ(θ, p, integrator)
.θ₀
: The estimates estimates for the MLEs.integrator
: The integrator for the differential equation problem. See alsoconstruct_integrator
.
Keyword Arguments
syms=eachindex(θ₀)
: Names for each parameter.data=SciMLBase.NullParameters()
: The parameterp
in the log-likelihood function.f_kwargs=nothing
: Keyword arguments, passed as aNamedTuple
, for theOptimizationFunction
.prob_kwargs=nothing
: Keyword arguments, passed as aNamedTuple
, for theOptimizationProblem
.
Outputs
Returns the LikelihoodProblem
problem object.
LikelihoodSolution
ProfileLikelihood.AbstractLikelihoodSolution
— Typeabstract type AbstractLikelihoodSolution{N, P}
Type representing the solution to a likelihood problem, where N
is the number of parameters and P
is the type of the likelihood problem.
ProfileLikelihood.LikelihoodSolution
— Typestruct LikelihoodSolution{Θ,P,M,R,A} <: AbstractLikelihoodSolution
Struct for a solution to a LikelihoodProblem
.
Fields
mle::Θ
: The MLEs.problem::P
: TheLikelihoodProblem
.optimiser::A
: The algorithm used for solving the optimisation problem.maximum::M
: The maximum likelihood.retcode::R
: TheSciMLBase.ReturnCode
.
Constructors
LikelihoodSolution(sol::SciMLBase.OptimizationSolution, prob::AbstractLikelihoodProblem; alg=sol.alg)
Constructs the likelihood solution from a solution to an OptimizationProblem
with a given LikelihoodProblem
.
ProfileLikelihood.mle
— Functionmle(prob::LikelihoodProblem, alg, args...; kwargs...)
mle(prob::LikelihoodProblem, alg::Tuple, args...; kwargs...)
Given the likelihood problem prob
and an optimiser alg
, finds the MLEs and returns a LikelihoodSolution
object. Extra arguments and keyword arguments for solve
can be passed through args...
and kwargs...
.
If alg
is a Tuple
, then the problem is re-optimised after each algorithm with the next element in alg, starting from alg[1]
, with initial estimate coming from the solution with the previous algorithm (starting with get_initial_estimate(prob)
).
ProfileLikelihoodSolution
ProfileLikelihood.ProfileLikelihoodSolution
— TypeProfileLikelihoodSolution{I,V,LP,LS,Spl,CT,CF,OM}
Struct for the normalised profile log-likelihood. See profile
for a constructor.
Fields
parameter_values::Dict{I, V}
This is a dictionary such that parameter_values[i]
gives the parameter values used for the normalised profile log-likelihood of the i
th variable.
profile_values::Dict{I, V}
This is a dictionary such that profile_values[i]
gives the values of the normalised profile log-likelihood function at the corresponding values in θ[i]
.
likelihood_problem::LP
The original LikelihoodProblem
.
likelihood_solution::LS
The solution to the full problem.
splines::Dict{I, Spl}
This is a dictionary such that splines[i]
is a spline through the data (parameter_values[i], profile_values[i])
. This spline can be evaluated at a point ψ
for the i
th variable by calling an instance of the struct with arguments (ψ, i)
. See also spline_profile
.
confidence_intervals::Dict{I,ConfidenceInterval{CT,CF}}
This is a dictonary such that confidence_intervals[i]
is a confidence interval for the i
th parameter.
other_mles::OM
This is a dictionary such that other_mles[i]
gives the vector for the MLEs of the other parameters not being profiled, for each datum.
Spline evaluation
This struct is callable. We define the method
(prof::ProfileLikelihoodSolution)(θ, i)
that evaluates the spline through the i
th profile at the point θ
.
ProfileLikelihood.ConfidenceInterval
— Typestruct ConfidenceInterval{T, F}
Struct representing a confidence interval.
Fields
lower::T
The lower bound of the confidence interval.
upper::T
The upper bound of the confidence interval.
level::F
The level of the confidence interval.
ProfileLikelihood.profile
— Functionprofile(prob::LikelihoodProblem, sol::LikelihoodSolution, n=1:number_of_parameters(prob);
alg=get_optimiser(sol),
conf_level::F=0.95,
confidence_interval_method=:spline,
threshold=get_chisq_threshold(conf_level),
resolution=200,
param_ranges=construct_profile_ranges(sol, get_lower_bounds(prob), get_upper_bounds(prob), resolution),
min_steps=10,
normalise::Bool=true,
spline_alg=FritschCarlsonMonotonicInterpolation,
extrap=Line,
parallel=false,
next_initial_estimate_method = :prev,
kwargs...)
Computes profile likelihoods for the parameters from a likelihood problem prob
with MLEs sol
.
See also replace_profile!
which allows you to re-profile a parameter in case you are not satisfied with the results. For plotting, see the plot_profiles
function (requires that you have loaded a backend of Makie.jl).
Arguments
prob::LikelihoodProblem
: TheLikelihoodProblem
.sol::LikelihoodSolution
: TheLikelihoodSolution
. See alsomle
.n=1:number_of_parameters(prob)
: The parameter indices to compute the profile likelihoods for.
Keyword Arguments
alg=get_optimiser(sol)
: The optimiser to use for solving each optimisation problem.conf_level::F=0.95
: The level to use for theConfidenceInterval
s.confidence_interval_method=:spline
: The method to use for computing the confidence intervals. See alsoget_confidence_intervals!
. The default:spline
uses rootfinding on the spline through the data, defining a continuous function, while the alternative:extrema
simply takes the extrema of the values that exceed the threshold.threshold=get_chisq_threshold(conf_level)
: The threshold to use for defining the confidence intervals.resolution=200
: The number of points to use for evaluating the profile likelihood in each direction starting from the MLE (giving a total of2resolution
points). -resolution=200
: The number of points to use for defininggrids
below, giving the number of points to the left and right of each interest parameter. This can also be a vector, e.g.resolution = [20, 50, 60]
will use20
points for the first parameter,50
for the second, and60
for the third.param_ranges=construct_profile_ranges(sol, get_lower_bounds(prob), get_upper_bounds(prob), resolution)
: The ranges to use for each parameter.min_steps=10
: The minimum number of steps to allow for the profile in each direction. If fewer than this number of steps are used before reaching the threshold, then the algorithm restarts and computes the profile likelihood a numbermin_steps
of points in that direction. See alsomin_steps_fallback
.min_steps_fallback=:replace
: Method to use for updating the profile when it does not reach the minimum number of steps,min_steps
. See alsoreach_min_steps!
. If:replace
, then the profile is completely replaced and we usemin_steps
equally spaced points to replace it. If:refine
, we just fill in some of the space in the grid so that amin_steps
number of points are reached. Note that this latter option will mean that the spacing is no longer constant between parameter values. You can use:refine_parallel
to apply:refine
in parallel.normalise::Bool=true
: Whether to optimise the normalised profile log-likelihood or not.spline_alg=FritschCarlsonMonotonicInterpolation
: The interpolation algorithm to use for computing a spline from the profile data. See Interpolations.jl.extrap=Line
: The extrapolation algorithm to use for computing a spline from the profile data. See Interpolations.jl.parallel=false
: Whether to use multithreading. Iftrue
, will use multithreading so that multiple parameters are profiled at once, and the steps to the left and right are done at the same time.next_initial_estimate_method = :prev
: Method for selecting the next initial estimate when stepping forward when profiling.:prev
simply uses the previous solution, but you can also use:interp
to use linear interpolation. See alsoset_next_initial_estimate!
.kwargs...
: Extra keyword arguments to pass intosolve
for solving theOptimizationProblem
. See also the docs from Optimization.jl.
Output
Returns a ProfileLikelihoodSolution
.
ProfileLikelihood.replace_profile!
— Functionreplace_profile!(prof::ProfileLikelihoodSolution, n);
alg=get_optimiser(prof.likelihood_solution),
conf_level::F=0.95,
confidence_interval_method=:spline,
threshold=get_chisq_threshold(conf_level),
resolution=200,
param_ranges=construct_profile_ranges(prof.likelihood_solution, get_lower_bounds(prof.likelihood_problem), get_upper_bounds(prof.likelihood_problem), resolution),
min_steps=10,
min_steps_fallback=:replace,
normalise::Bool=true,
spline_alg=FritschCarlsonMonotonicInterpolation,
extrap=Line,
parallel=false,
next_initial_estimate_method=:prev,
kwargs...) where {F}
Given an existing prof::ProfileLikelihoodSolution
, replaces the profile results for the parameters in n
by re-profiling. The keyword arguments are the same as for profile
.
ProfileLikelihood.refine_profile!
— Functionrefine_profile!(prof::ProfileLikelihoodSolution, n;
alg=get_optimiser(prof.likelihood_solution),
conf_level::F=0.95,
confidence_interval_method=:spline,
threshold=get_chisq_threshold(conf_level),
target_number=10,
normalise::Bool=true,
spline_alg=FritschCarlsonMonotonicInterpolation,
extrap=Line,
parallel=false,
kwargs...) where {F}
Given an existing prof::ProfileLikelihoodSolution
, refines the profile results for the parameters in n
by adding more points. The keyword arguments are the same as for profile
. target_number
is the total number of points that should be included in the end (not how many more are added).
ProfileLikelihood.set_next_initial_estimate!
— Methodset_next_initial_estimate!(sub_cache, param_vals, other_mles, prob, θₙ; next_initial_estimate_method=Val(:prev))
Method for selecting the next initial estimate for the optimisers. sub_cache
is the cache vector for placing the initial estimate into, param_vals
is the current list of parameter values for the interest parameter, and other_mles
is the corresponding list of previous optimisers. prob
is the OptimizationProblem
. The value θₙ
is the next value of the interest parameter.
The available methods are:
next_initial_estimate_method = Val(:prev)
: If this is selected, simply useother_mles[end]
, i.e. the previous optimiser.next_initial_estimate_method = Val(:interp)
: If this is selected, the next optimiser is determined via linear interpolation using the data(param_vals[end-1], other_mles[end-1]), (param_vals[end], other_mles[end])
. If the new approximation is outside of the parameter bounds, falls back tonext_initial_estimate_method = :prev
.
ProfileLikelihood.get_confidence_intervals!
— Functionget_confidence_intervals!(confidence_intervals, method, n, param_vals, profile_vals, threshold, spline_alg, extrap, mles, conf_level)
Method for computing the confidence intervals.
Arguments
confidence_intervals
: The dictionary storing the confidence intervals.method
: The method to use for computing the confidence interval. The available methods are:method = Val(:spline)
: Fits a spline to(param_vals, profile_vals)
and finds where the continuous spline equalsthreshold
.method = Val(:extrema)
: Takes the first and last values inparam_vals
whose corresponding value inprofile_vals
exceedsthreshold
.
n
: The parameter being profiled.param_vals
: The parameter values.profile_vals
: The profile values.threshold
: The threshold for the confidence interval.spline_alg
: The algorithm to use for fitting a spline.extrap
: The extrapolation algorithm used for the spline.mles
: The MLEs.conf_level
: The confidence level for the confidence interval.
Outputs
There are no outputs - confidence_intervals[n]
gets the ConfidenceInterval
put into it.
ProfileLikelihood.reach_min_steps!
— Functionreach_min_steps!(param_vals, profile_vals, other_mles, param_range,
restricted_prob, n, cache, alg, sub_cache, ℓmax, normalise,
threshold, min_steps, mles; min_steps_fallback=Val(:replace), next_initial_estimate_method=Val(:interp), kwargs...)
Updates the results from the side of a profile likelihood (e.g. left or right side, see find_endpoint!
) to meet the minimum number of steps min_steps
.
Arguments
param_vals
: The parameter values.profile_vals
: The profile values.other_mles
: The other MLEs, i.e. the optimised parameters for the corresponding fixed parameter values inparam_vals
.param_range
: The vector of parameter values.restricted_prob
: The optimisation problem, restricted to then
th parameter.n
: The parameter being profiled.cache
: A cache for the complete parameter vector.alg
: The algorithm used for optimising.sub_cache
: A cache for the parameter vector excluding then
th parameter.ℓmax
: The maximum likelihood.normalise
: Whether the optimisation problem is normalised.threshold
: The threshold for the confidence interval.min_steps
: The minimum number of steps to reach.mles
: The MLEs.
Keyword Arguments
min_steps_fallback=Val(:interp)
: The method used for reaching the minimum number of steps. The available methods are:min_steps_fallback = Val(:replace)
: This method completely replaces the profile, defining a grid from the MLE to the computed endpoint withmin_steps
points. No information is re-used.min_steps_fallback = Val(:refine)
: This method just adds more points to the profile, filling in enough points so that the total number of points ismin_steps
. The initial estimates in this case come from a spline fromother_mles
.min_steps_fallback = Val(:parallel_refine)
: This applies the method above, except in parallel.
next_initial_estimate_method=Val(:replace)
: The method used for obtaining initial estimates. See alsoset_next_initial_estimate!
.
BivariateProfileLikelihoodSolution
ProfileLikelihood.BivariateProfileLikelihoodSolution
— TypeBivariateProfileLikelihoodSolution{I,V,LP,LS,Spl,CT,CF,OM}
Struct for the normalised bivariate profile log-likelihood. See bivariate_profile
for a constructor.
Arguments
parameter_values::Dict{I, G}
Maps the tuple (i, j)
to the grid values used for this parameter pair. The result is a Tuple
, with the first element the grid for the i
th parameter, and the second element the grid for the j
th parameter. The grids are given as OffsetVector
s, with the 0
th index the MLE, negative indices to the left of the MLE, and positive indices to the right of the MLE.
profile_values::Dict{I, V}
Maps the tuple (i, j)
to the matrix used for this parameter pair. The result is a OffsetMatrix
, with the (k, ℓ)
entry the profile at (parameter_values[(i, j)][1][k], parameter_values[(i, j)][2][k])
, and particularly the (0, 0)
entry is the profile at the MLEs.
likelihood_problem::LP
The original likelihood problem.
likelihood_solution::LS
The original likelihood solution.
interpolants::Dict{I,Spl}
Maps the tuple (i, j)
to the interpolant for that parameter pair's profile. This interpolant also uses linear extrapolation.
confidence_regions::Dict{I,ConfidenceRegion{CT,CF}}
Maps the tuple (i, j)
to the confidence region for that parameter pair's confidence region. See also ConfidenceRegion
.
other_mles::OM
Maps the tuple (i, j)
to an OffsetMatrix
storing the solutions for the nuisance parameters at the corresponding grid values.
Interpolant evaluation
This struct is callable. We define the method
(prof::BivariateProfileLikelihoodSolution)(θ, ψ, i, j)
that evaluates the interpolant through the (i, j)
th profile at the point (θ, ψ)
.
ProfileLikelihood.ConfidenceRegion
— Typestruct ConfidenceRegion{T, F}
Struct representing a confidence region.
Fields
x::T
The x
-coordinates for the region's boundary.
y::T
The y
-coordinates for the region's boundary.
level::F
The level of the confidence region.
ProfileLikelihood.bivariate_profile
— Functionbivariate_profile(prob::LikelihoodProblem, sol::LikelihoodSolution, n::NTuple{M,NTuple{2,Int}};
alg=get_optimiser(sol),
conf_level::F=0.95,
confidence_region_method=Val(:contour),
threshold=get_chisq_threshold(conf_level, 2),
resolution=200,
grids=construct_profile_grids(n, sol, get_lower_bounds(prob), get_upper_bounds(prob), resolution),
min_layers=10,
outer_layers=0,
normalise=Val(true),
parallel=Val(false),
next_initial_estimate_method=Val(:nearest),
kwargs...) where {M,F}
Computes bivariates profile likelihoods for the parameters from a likelihood problem prob
with MLEs sol
. You can also call this function using Symbols
, e.g. if get_syms(prob) = [:λ, :K, :u₀]
, then calling bivariate_profile(prob, sol, ((:λ, :K), (:K, u₀)))
is the same as calling bivariate_profile(prob, sol, ((1, 2), (2, 3)))
(the integer coordinate representation is still used in the solution, though).
For plotting, see the plot_profiles
function (requires that you have loaded a backend of Makie.jl).
Arguments
prob::LikelihoodProblem
: TheLikelihoodProblem
.sol::LikelihoodSolution
: TheLikelihoodSolution
. See alsomle
.n::NTuple{M,NTuple{2,Int}}
: The parameter indices to compute the profile likelihoods for. These should be tuples of indices, e.g.n = ((1, 2),)
will compute the bivariate profile between the parameters1
and2
.
Keyword Arguments
alg=get_optimiser(sol)
: The optimiser to use for solving each optimisation problem.conf_level::F=0.95
: The level to use for theConfidenceRegion
s.confidence_region_method=:contour
: The method to use for computing the confidence regions. See alsoget_confidence_regions!
. The default,:contour
, using Contour.jl to compute the boundary of the confidence region. An alternative option is:delaunay
, which uses DelaunayTriangulation.jl and triangulation contouring to find the boundary. This latter option is only available if you have already doneusing DelaunayTriangulation
.threshold=get_chisq_threshold(conf_level, 2)
: The threshold to use for defining the confidence regions.resolution=200
: The number of points to use for defininggrids
below, giving the number of points to the left and right of each interest parameter. This can also be a vector, e.g.resolution = [20, 50, 60]
will use20
points for the first parameter,50
for the second, and60
for the third. When defining the grid between pairs of values, the maximum of the two resolutions is used (thus defining a square grid).grids=construct_profile_grids(n, sol, get_lower_bounds(prob), get_upper_bounds(prob), resolution)
: The grids to use for each parameter pair.min_layers=10
: The minimum number of layers to allow for the profile away from the MLE. If fewer than this number of layers are used before reaching the threshold, then the algorithm restarts and computes the profile likelihood a numbermin_steps
of points in that direction.outer_layers=0
: The number of layers to go out away from the bounding box of the confidence region.normalise=true
: Whether to optimise the normalised profile log-likelihood or not.parallel=false
: Whether to use multithreading. Iftrue
, will use multithreading so that multiple parameters are profiled at once, and the work done evaluating the solution at each node in a layer is distributed across each thread.next_initial_estimate_method = :nearest
: Method for selecting the next initial estimate when stepping onto the next layer when profiling.:nearest
simply uses the solution at the nearest node from the previous layer, but you can also use:mle
to reuse the MLE or:interp
to use linear interpolation. See alsoset_next_initial_estimate!
.kwargs...
: Extra keyword arguments to pass intosolve
for solving theOptimizationProblem
. See also the docs from Optimization.jl.
Output
Returns a BivariateProfileLikelihoodSolution
.
ProfileLikelihood.set_next_initial_estimate!
— Methodset_next_initial_estimate!(sub_cache, other_mles, I::CartesianIndex, fixed_vals, grid, layer, prob, next_initial_estimate_method::Val{M}) where {M}
Method for selecting the next initial estimate for the optimisers.
Arguments
sub_cache
: Cache for the next initial estimate.other_mles
: Solutions to the optimisation problems found so far.I::CartesianIndex
: The coordinate of the node currently being considered.fixed_vals
: The current values for the parameters of interest.grid
: The grid for the parameters of interest.layer
: The current layer.prob
: The restricted optimisation problem.next_initial_estimate_method::Val{M}
: The method to use.
The methods currently available for next_initial_estimate_method
are:
next_initial_estimate_method = Val(:mle)
: Simply setssub_cache
to be the MLE.next_initial_estimate_method = Val(:nearest)
: Setssub_cache
to beother_mles[J]
, whereJ
is the nearest node toI
in the previous layer.next_initial_estimate_method = Val(:interp)
: Uses linear interpolation from all the previous layers to extrapolate and compute a newsub_cache
.
Outputs
There are no outputs.
Prediction intervals
ProfileLikelihood.get_prediction_intervals
— Functionget_prediction_intervals(q, prof::(Bivariate)ProfileLikelihoodSolution, data;
q_prototype=isinplace(q, 3) ? nothing : build_q_prototype(q, prof, data), resolution=250)
Obtain prediction intervals for the output of the prediction function q
, assuming q
returns (or operates in-place on) a vector.
Arguments
q
: The prediction function, taking either the form(θ, data)
or(cache, θ, data)
. The former version is an out-of-place version, returning the full vector, while the latter version is an in-place version, with the output being placed intocache
. The argumentθ
is the same as the parameters used in the likelihood problem (fromprof
), and thedata
argument is the samedata
as in this function.prof::(Bivariate)ProfileLikelihoodSolution
: The profile likelihood results.data
: The argumentdata
inq
.
Keyword Arguments
q_prototype=isinplace(q, 3) ? nothing : build_q_prototype(q, prof, data)
: A prototype for the result ofq
. If you are using theq(θ, data)
version ofq
, this can be inferred frombuild_q_prototype
, but if you are using the in-place version then abuild_q_prototype
is needed. For example, ifq
returns a vector with eltypeFloat64
and has length 27,q_prototype
could bezeros(27)
.resolution::Integer=250
: The amount of curves to evaluate for each parameter. This will be the same for each parameter. Ifprof isa BivariateProfileLikelihoodSolution
, thenresolution^2
points are defined inside a bounding box for the confidence region, and then we throw away all points outside of the actual confidence region.parallel=false
: Whether to use multithreading. Multithreading is used when buildingq_vals
below.
Outputs
Four values are returned. In order:
individual_intervals
: Prediction intervals for the output ofq
, relative to each parameter.union_intervals
: The union of the individual prediction intervals fromindividual_intervals
.q_vals
: Values ofq
at each parameter considered. The output is aDict
, where the parameter index is mapped to a matrix where each column is an output fromq
, with thej
th column corresponding to the parameter value atparam_ranges[j]
.param_ranges
: Parameter values used for each prediction interval.
Plotting
ProfileLikelihood.plot_profiles
— Functionplot_profiles(prof::ProfileLikelihoodSolution, vars = profiled_parameters(prof);
ncol=nothing,
nrow=nothing,
true_vals=Dict(vars .=> nothing),
spline=true,
show_mles=true,
shade_ci=true,
fig_kwargs=nothing,
axis_kwargs=nothing,
show_points=false,
markersize=9,
latex_names = Dict(vars .=> [L" heta_{i}" for i in SciMLBase.sym_to_index.(vars, Ref(prof))]))
Plot results from a profile likelihood solution prof
. To use this function you you need to have done using CairoMakie
(or any other Makie backend).
Arguments
prof::ProfileLikelihoodSolution
: The profile likelihood solution fromprofile
.vars = profiled_parameters(prof)
: The parameters to plot.
Keyword Arguments
ncol=nothing
: The number of columns to use. Ifnothing
, chosen automatically viachoose_grid_layout
.nrow=nothing
: The number of rows to use. Ifnothing
, chosen automatically viachoose_grid_layout
true_vals=Dict(vars .=> nothing)
: A dictionary mapping parameter indices to their true values, if they exist. Ifnothing
, nothing is plotted, otherwise a black line is plotted at the true value for the profile.spline=true
: Whether the curve plotted should come from a spline through the results, or if the data itself should be plotted.show_mles=true
: Whether to put a red line at the MLEs.shade_ci=true
: Whether to shade the area under the profile between the confidence interval.fig_kwargs=nothing
: Extra keyword arguments forFigure
(see the Makie docs).axis_kwargs=nothing
: Extra keyword arguments forAxis
(see the Makie docs).show_points=false
: Whether to show the profile data.markersize=9
: The marker size used forshow_points
.latex_names = Dict(vars .=> [L" heta_{i}" for i in SciMLBase.sym_to_index.(vars, Ref(prof))]))
: LaTeX names to use for the parameters. Defaults toθᵢ
, wherei
is the index of the parameter.
Output
The Figure()
is returned.
plot_profiles(prof::BivariateProfileLikelihoodSolution, vars = profiled_parameters(prof);
ncol=nothing,
nrow=nothing,
true_vals=Dict(1:number_of_parameters(get_likelihood_problem(prof)) .=> nothing),
show_mles=true,
fig_kwargs=nothing,
axis_kwargs=nothing,
interpolation=false,
smooth_confidence_boundary=false,
close_contour=true,
latex_names=Dict(1:number_of_parameters(get_likelihood_problem(prof)) .=> get_syms(prof)))
Plot results from a bivariate profile likelihood solution prof
. To use this function you you need to have done using CairoMakie
(or any other Makie backend).
Arguments
prof::ProfileLikelihoodSolution
: The profile likelihood solution fromprofile
.vars = profiled_parameters(prof)
: The parameters to plot.
Keyword Arguments
ncol=nothing
: The number of columns to use. Ifnothing
, chosen automatically viachoose_grid_layout
.nrow=nothing
: The number of rows to use. Ifnothing
, chosen automatically viachoose_grid_layout
true_vals=Dict(1:number_of_parameters(get_likelihood_problem(prof)) .=> nothing)
: A dictionary mapping parameter indices to their true values, if they exist. Ifnothing
, nothing is plotted, otherwise a black dot is plotted at the true value on the bivariate profile's plot.show_mles=true
: Whether to put a red dot at the MLEs.fig_kwargs=nothing
: Extra keyword arguments forFigure
(see the Makie docs).axis_kwargs=nothing
: Extra keyword arguments forAxis
(see the Makie docs).interpolation=false
: Whether to plot the profile using the interpolant (true
), or to use the data fromprof
directly (false
).smooth_confidence_boundary=false
: Whether to smooth the confidence region boundary when plotting (true
) or not (false
). The smoothing is done with a spline.close_contour=true
: Whether to connect the last part of the confidence region boundary to the beginning (true
) or not (false
).latex_names=Dict(1:number_of_parameters(get_likelihood_problem(prof)) .=> get_syms(prof)))
: LaTeX names to use for the parameters. Defaults to thesyms
names.xlim_tuples=nothing
:xlims
to use for each plot.nothing
if thexlims
should be set automatically.ylim_tuples=nothing
:ylims
to use for each plot.nothing
if theylims
should be set automatically.
Output
The Figure()
is returned.
GridSearch
Grid definitions
ProfileLikelihood.AbstractGrid
— Typeabstract type AbstractGrid{N,B,T}
Type representing a grid, where N
is the number of parameters, B
is the type for the bounds, and T
is the number type.
ProfileLikelihood.RegularGrid
— Typestruct RegularGrid{N,B,R,S,T} <: AbstractGrid{N,B,T}
Struct for a grid in which each parameter is regularly spaced.
Fields
lower_bounds::B
: Lower bounds for each parameter.upper_bounds::B
: Upper bounds for each parameter.resolution::R
: Number of grid points for each parameter. IfR <: Number
, then the same number of grid points is used for each parameter.step_sizes::S
: Grid spacing for each parameter.
Constructor
You can construct a RegularGrid
using RegularGrid(lower_bounds, upper_bounds, resolution)
.
ProfileLikelihood.FusedRegularGrid
— Typestruct FusedRegularGrid{N,B,R,S,T,C,OR} <: AbstractGrid{N,B,T}
Struct representing the fusing of two grids.
Fields
positive_grid::RegularGrid{N,B,R,S,T}
This is the first part of the grid, indexed into by positive integers.
negative_grid::RegularGrid{N,B,R,S,T}
This is the second part of the grid, indexed into by negative integers.
centre::C
The two grids meet at a common centre, and this is that centre
.
resolutions::R
This is the vector of resolutions provided (e.g. if store_original_resolutions=true
in the constructor below), or the transformed version from get_resolution_tuples
.
Constructor
You can construct a FusedRegularGrid
using the method
FusedRegularGrid(lower_bounds::B, upper_bounds::B, centre::C, resolutions::R; store_original_resolutions=false) where {B,R,C}
For example, the following code creates fused
as the fusion of grid_1
and grid_2
:
lb = [2.0, 3.0, 1.0, 5.0]
ub = [15.0, 13.0, 27.0, 10.0]
centre = [7.3, 8.3, 2.3, 7.5]
grid_1 = RegularGrid(centre .+ (ub .- centre) / 173, ub, 173)
grid_2 = RegularGrid(centre .- (centre .- lb) / 173, lb, 173)
fused = ProfileLikelihood.FusedRegularGrid(lb, ub, centre, 173)
There are 173
points to the left and right of centre
in this case. To use a varying number of points, use e.g.
lb = [2.0, 3.0, 1.0, 5.0, 4.0]
ub = [15.0, 13.0, 27.0, 10.0, 13.0]
centre = [7.3, 8.3, 2.3, 7.5, 10.0]
res = [(2, 11), (23, 25), (19, 21), (50, 51), (17, 99)]
grid_1 = RegularGrid(centre .+ (ub .- centre) ./ [2, 23, 19, 50, 17], ub, [2, 23, 19, 50, 17])
grid_2 = RegularGrid(centre .- (centre .- lb) ./ [11, 25, 21, 51, 99], lb, [11, 25, 21, 51, 99])
fused = ProfileLikelihood.FusedRegularGrid(lb, ub, centre, res) # fused grid_1 and grid_2
ProfileLikelihood.IrregularGrid
— Typestruct IrregularGrid{N,B,R,S,T} <: AbstractGrid{N,B,T}
Struct for an irregular grid of parameters.
Fields
lower_bounds::B
: Lower bounds for each parameter.upper_bounds::B
: Upper bounds for each parameter.grid::G
: The set of parameter values, e.g. a matrix where each column is the parameter vector.
Constructor
You can construct a IrregularGrid
using IrregularGrid(lower_bounds, upper_bounds, grid)
.
Performing a grid search
ProfileLikelihood.GridSearch
— Typestruct GridSearch{F,G}
Struct for a GridSearch
.
Fields
f::F
: The function to optimise, of the formf(x, p)
.p::P
: The argumentsp
in the functionf
.grid::G
: The grid, whereG<:AbstractGrid
. See alsogrid_search
.
ProfileLikelihood.grid_search
— Functiongrid_search(prob; save_vals=Val(false), minimise:=Val(false), parallel=Val(false))
Performs a grid search for the given grid search problem prob
.
Arguments
prob::GridSearch{F, G}
: The grid search problem.
Keyword Arguments
save_vals:=Val(false)
: Whether to return a array with the function values.minimise:=Val(false)
: Whether to minimise or to maximise the function.parallel:=Val(false)
: Whether to run the grid search with multithreading.
Outputs
f_opt
: The optimal objective value.x_argopt
: The parameter that gavef_opt
.f_res
: Ifsave_vals==Val(true)
, then this is the array of function values.
grid_search(f, grid::AbstractGrid; save_vals=Val(false), minimise=Val(false), parallel=Val(false))
For a given grid
and function f
, performs a grid search.
Arguments
f
: The function to optimise.grid::AbstractGrid
: The grid to use for optimising.
Keyword Arguments
save_vals=Val(false)
: Whether to return a array with the function values.minimise=Val(false)
: Whether to minimise or to maximise the function.parallel=Val(false)
: Whether to run the grid search with multithreading.
grid_search(prob::LikelihoodProblem, grid::AbstractGrid, parallel=Val(false); save_vals=Val(false))
Given a grid
and a likelihood problem prob
, maximises it over the grid using a grid search. If save_vals==Val(true)
, then the likelihood function values at each gridpoint are returned. Set parallel=Val(true)
if you want multithreading.