API references

The package exports the following functions for parameters identifiability analysis, confidence intervals evaluation and results visualization.

LikelihoodProfiler.get_endpointFunction
function get_endpoint(
    theta_init::Vector{Float64},
    scan_func::Function,
    loss_func::Function,
    method::Symbol,
    direction::Symbol = :right;

    loss_crit::Float64 = 0.0,
    scale::Vector{Symbol} = fill(:direct, length(theta_init)),
    theta_bounds::Vector{Tuple{Float64,Float64}} = unscaling.(
        fill((-Inf, Inf), length(theta_init)),
        scale
        ),
    scan_bound::Float64 = unscaling(
        (direction==:left) ? -9.0 : 9.0,
        scale[theta_num]
        ),
    scan_tol::Float64 = 1e-3,
    loss_tol::Float64 = 0.,
    local_alg::Symbol = :LN_NELDERMEAD,
    kwargs...
    )

Calculates confidence interval's right or left endpoints for a function of parameters scan_func.

Return

EndPoint object storing confidence interval's endpoint and intermediate profile points.

Arguments

  • theta_init: starting values of parameter vector $\theta$. The starting values should not necessary be the optimum values of loss_func but loss_func(theta_init) should be lower than loss_crit.
  • scan_func: function of parameters.
  • loss_func: loss function $\Lambda\left(\theta\right)$ for profile likelihood-based (PL) identification. Usually we use log-likelihood for PL analysis: $\Lambda( \theta ) = - 2 ln\left( L(\theta) \right)$.
  • method: computational method to estimate confidence interval's endpoint. Currently the only supported method is: :CICO_ONE_PASS.
  • direction: :right or :left endpoint to estimate.

Keyword arguments

  • loss_crit: critical level of loss function. Confidence interval's endpoint value is the intersection point of profile likelihood and loss_crit level.
  • scale: vector of scale transformations for each parameters' component. Possible values: :direct (:lin), :log, :logit. This option can speed up the optimization, especially for wide theta_bounds. The default value is :direct (no transformation) for all parameters.
  • theta_bounds: vector of tuple (lower_bound, upper_bound) for each parameter. Bounds define the ranges for possible parameter values. Default bounds are (-Inf,Inf).
  • scan_bound: value which states the area of confidence point analysis.
  • scan_tol: Absolute tolerance for theta_num parameter used as termination criterion.
  • loss_tol: Absolute tolerance controlling loss_func closeness to loss_crit (termination criterion). Currently doesn't work for :CICO_ONE_PASS method because of limitation in LN_AUGLAG interface.
  • local_alg: algorithm of optimization. Derivative-free and gradient-based algorithms form NLopt package.
  • max_iter : maximal number of fitter iterations. If reaches the result status will be :MAX_ITER_STOP.
  • scan_grad : For gradient optimization methods it is necessary to set how the gradient of scan_func should be calculated.
    • :EMPTY (default) no gradient is set. It works only for gradient-free methods.
    • :AUTODIFF means autodifferentiation from ForwardDiff package is used.
    • :FINITE means finite difference method from Calculus is used.
    • function(x::Vector{Float64}) which returns gradient vector.
  • loss_grad : For gradient optimization methods it is necessary to set how the gradient of loss_func should be calculated.
    • :EMPTY (default) no gradient is set. It works only for gradient-free methods.
    • :AUTODIFF means autodifferentiation from ForwardDiff package is used.
    • :FINITE means finite difference method from Calculus is used.
    • function(x::Vector{Float64}) which returns gradient vector.
  • silent : Boolean argument declaring whether we display the optimization progress. Default is false

See also get_interval

source
LikelihoodProfiler.get_endpointFunction
function get_endpoint(
    theta_init::Vector{Float64},
    theta_num::Int,
    loss_func::Function,
    method::Symbol,
    direction::Symbol = :right;

    loss_crit::Float64 = 0.0,
    scale::Vector{Symbol} = fill(:direct, length(theta_init)),
    theta_bounds::Vector{Tuple{Float64,Float64}} = unscaling.(
        fill((-Inf, Inf), length(theta_init)),
        scale
        ),
    scan_bound::Float64 = unscaling(
        (direction==:left) ? -9.0 : 9.0,
        scale[theta_num]
        ),
    scan_tol::Float64 = 1e-3,
    loss_tol::Float64 = 0.,
    local_alg::Symbol = :LN_NELDERMEAD,
    max_iter::Int = 10^5,
    loss_grad::Union{Function, Symbol} = :EMPTY,
    silent::Bool = false
)

Calculates confidence interval's right or left endpoints for a given parameter theta_num.

Return

EndPoint object storing confidence interval's endpoint and intermediate profile points.

Arguments

  • theta_init: starting values of parameter vector $\theta$. The starting values should not necessary be the optimum values of loss_func but loss_func(theta_init) should be lower than loss_crit.
  • theta_num: index of vector component for identification: theta_init(theta_num).
  • loss_func: loss function $\Lambda\left(\theta\right)$ for profile likelihood-based (PL) identification. Usually we use log-likelihood for PL analysis: $\Lambda( \theta ) = - 2 ln\left( L(\theta) \right)$.
  • method: computational method to estimate confidence interval's endpoint. Currently the following methods are implemented: :CICO_ONE_PASS, :LIN_EXTRAPOL, :QUADR_EXTRAPOL.
  • direction: :right or :left endpoint to estimate.

Keyword arguments

  • loss_crit : critical level of loss function. Confidence interval's endpoint value is the intersection point of profile likelihood and loss_crit level.
  • scale : vector of scale transformations for each parameters' component. Possible values: :direct (:lin), :log, :logit. This option can speed up the optimization, especially for wide theta_bounds. The default value is :direct (no transformation) for all parameters.
  • theta_bounds : vector of tuple (lower_bound, upper_bound) for each parameter. Bounds define the ranges for possible parameter values. Default bounds are (-Inf,Inf).
  • scan_bound : value which states the area of confidence point analysis.
  • scan_tol : Absolute tolerance for theta_num parameter used as termination criterion.
  • loss_tol : Absolute tolerance controlling loss_func closeness to loss_crit (termination criterion). Currently doesn't work for :CICO_ONE_PASS method because of limitation in LN_AUGLAG interface.
  • local_alg : algorithm of optimization. Derivative-free and gradient-based algorithms form NLopt package.
  • max_iter : maximal number of fitter iterations. If reaches the result status will be :MAX_ITER_STOP.
  • loss_grad : For gradient optimization methods it is necessary to set how the gradient of loss_func should be calculated. There are options:
    • :EMPTY (default) no gradient is set. It works only for gradient-free methods.
    • :AUTODIFF means autodifferentiation from ForwardDiff package is used.
    • :FINITE means finite difference method from Calculus is used.
    • It is also possible to set gradient function here function(x::Vector{Float64}) which returns gradient vector.
  • silent : Boolean argument declaring whether we display the optimization progress. Default is false

See also get_interval

source
LikelihoodProfiler.get_intervalMethod
get_interval(
    theta_init::Vector{Float64},
    scan_func::Function,
    loss_func::Function,
    method::Symbol;

    loss_crit::Float64 = 0.0,
    scale::Vector{Symbol} = fill(:direct, length(theta_init)),
    theta_bounds::Vector{Tuple{Float64,Float64}} = unscaling.(
        fill((-Inf, Inf), length(theta_init)),
        scale
        ),
    scan_bounds::Tuple{Float64,Float64} = unscaling.(
        (-9.0, 9.0),
        :direct
        ),
    scan_tol::Float64 = 1e-3,
    loss_tol::Float64 = 0.,
    local_alg::Symbol = :LN_NELDERMEAD,
    autodiff::Bool = true,
    kwargs...
    )

Computes confidence interval for function of parameters scan_func.

Return

ParamInterval structure storing all input data and estimated confidence interval.

Arguments

  • theta_init: starting values of parameter vector $\theta$. The starting values should not necessary be the optimum values of loss_func but loss_func(theta_init) should be lower than loss_crit.
  • scan_func: scan function of parameters.
  • loss_func: loss function $\Lambda\left(\theta\right)$ the profile of which is analyzed. Usually we use log-likelihood for profile analysis in form $\Lambda( \theta ) = - 2 ln\left( L(\theta) \right)$.
  • method: computational method to estimate confidence interval's endpoint. Currently supports only :CICO_ONE_PASS method.

Keyword arguments

  • scale: vector of scale transformations for each parameters' component. Possible values: :direct (:lin), :log, :logit. This option can speed up the optimization, especially for wide theta_bounds. The default value is :direct (no transformation) for all parameters.
  • scan_bounds: scan bounds tuple for scan_func values. Default is (1e-9, 1e9) .
  • kwargs...: the additional arguments passed to get_endpoint
source
LikelihoodProfiler.get_intervalMethod
function get_interval(
    theta_init::Vector{Float64},
    theta_num::Int,
    loss_func::Function,
    method::Symbol;

    scale::Vector{Symbol} = fill(:direct, length(theta_init)),
    scan_bounds::Tuple{Float64,Float64} = unscaling.(
        (-9.0, 9.0),
        scale[theta_num]
        ),
    kwargs...
)

Computes confidence interval for single component theta_num of parameter vector.

Return

ParamInterval structure storing all input data and estimated confidence interval.

Arguments

  • theta_init: starting values of parameter vector $\theta$. The starting values should not necessary be the optimum values of loss_func but loss_func(theta_init) should be lower than loss_crit.
  • theta_num: index of vector component for identification: theta_init(theta_num).
  • loss_func: loss function $\Lambda\left(\theta\right)$ for profile likelihood-based (PL) identification. Usually we use log-likelihood for PL analysis: $\Lambda( \theta ) = - 2 ln\left( L(\theta) \right)$.
  • method: computational method to estimate confidence interval's endpoint. Currently the following methods are implemented: :CICO_ONE_PASS, :LIN_EXTRAPOL, :QUADR_EXTRAPOL.

Keyword arguments

  • scale: vector of scale transformations for each parameters' component. Possible values: :direct (:lin), :log, :logit. This option can speed up the optimization, especially for wide theta_bounds. The default value is :direct (no transformation) for all parameters.
  • scan_bounds: scan bounds tuple for theta_num parameter. Should be within the theta_bounds for theta_num parameter. Default is (-9.,9.) for :direct scales and (1e-9, 1e+9) for :log.
  • kwargs...: the additional arguments passed to get_endpoint
source
LikelihoodProfiler.get_optimalMethod
function get_optimal(
    theta_init::Vector{Float64}, # initial point of parameters
    loss_func::Function; # lambda(theta)

    scale::Vector{Symbol} = fill(:direct, length(theta_init)),
    theta_bounds::Vector{Tuple{Float64,Float64}} = unscaling.(
        fill((-Inf, Inf), length(theta_init)),
        scale
        ),
    scan_tol::Float64 = 1e-3,
    loss_tol::Float64 = 1e-3,
    local_alg::Symbol = :LN_NELDERMEAD,
    silent::Bool = false,
    max_iter::Int = 10^5,
    loss_grad::Union{Function, Symbol} = :EMPTY
)

Provides the optimization routine using the interface similar to get_endpoint. Currently it uses standard NLopt optimization but allows to use parameter scaling and autodiff method.

Return

ProfilePoint object storing optimizations results.valueandlossproperties are equal to the optimal (minimal)loss_func` value.

Arguments

  • theta_init: starting values of parameter vector $\theta$.
  • loss_func: loss function $\Lambda\left(\theta\right)$ for profile likelihood-based (PL) identification. Usually we use log-likelihood for PL analysis: $\Lambda( \theta ) = - 2 ln\left( L(\theta) \right)$.

Keyword arguments

  • scale: vector of scale transformations for each parameters' component. Possible values: :direct (:lin), :log, :logit. This option can speed up the optimization, especially for wide theta_bounds. The default value is :direct (no transformation) for all parameters.
  • theta_bounds: vector of tuple (lower_bound, upper_bound) for each parameter. Bounds define the ranges for possible parameter values. Default bounds are (-Inf,Inf).
  • scan_tol: Absolute tolerance for all component of theta used as termination criterion.
  • loss_tol: Absolute tolerance controlling loss_func.
  • local_alg: algorithm of optimization. Derivative-free and gradient-based algorithms form NLopt package.
  • max_iter : maximal number of optimizer iterations.
  • loss_grad : This option declares the method for calculating loss function gradient. This is required for gradient-based methods. The possible values - :EMPTY (default) not gradient is set. IT works only for gradient-free methods. - :AUTODIFF means autodifferentiation from ForwardDiff package is used. - :FINITE means finite difference method from Calculus is used. - It is also possible to set gradient function here function(x::Vector{Float64}) which returns gradient vector.
source
LikelihoodProfiler.profileMethod
function profile(
    theta_init::Vector{Float64},
    theta_num::Int,
    loss_func::Function;

    skip_optim::Bool = false,
    scale::Vector{Symbol} = fill(:direct, length(theta_init)),
    theta_bounds::Vector{Tuple{Float64,Float64}} = unscaling.(
        fill((-Inf, Inf), length(theta_init)),
        scale
        ),
    local_alg::Symbol = :LN_NELDERMEAD,
    ftol_abs::Float64 = 1e-3,
    maxeval::Int = 10^5,
    kwargs... # currently not used
)

Generates the profile function based on loss_func. Used internally in methods :LIN_EXTRAPOL, :QUADR_EXTRAPOL.

Return

Returns profile function for selected parameter.

Arguments

  • theta_init: starting values of parameter vector $\theta$.
  • theta_num: index of vector component for identification: theta_init(theta_num).
  • loss_func: loss function $\Lambda\left(\theta\right)$ for profile likelihood-based (PL) identification. Usually we use log-likelihood for PL analysis: $\Lambda( \theta ) = - 2 ln\left( L(\theta) \right)$.

Keyword arguments

  • skip_optim : set true if you need marginal profile, i.e. profile without optimization. Default is false.
  • scale :
  • theta_bounds: vector of tuple (lower_bound, upper_bound) for each parameter. Bounds define the ranges for possible parameter values. Default bounds are (-Inf,Inf).
  • local_alg: algorithm of optimization. Derivative-free and gradient-based algorithms form NLopt package.
  • ftol_abs : absolute tolerance criterion for profile function.
  • maxeval : maximal number of loss_func evaluations to estimate profile point.
source
LikelihoodProfiler.scalingFunction
scaling(x::Real, scale::Symbol = :direct)

Transforms values from specific scale to range [-Inf, Inf] based on option.

Return

Transformed value.

Arguments

  • x: input value.
  • scale: transformation type: :direct (:lin), :log, :logit.
source
LikelihoodProfiler.unscalingFunction
unscaling(x::Real, scale::Symbol = :direct)

Transforms values from [-Inf, Inf] to specific scale based on option. Inverse function for scaling.

Return

Transformed value.

Arguments

  • x: input value.
  • scale: transformation type: :direct (:lin), :log, :logit.
source
LikelihoodProfiler.update_profile_points!Method
update_profile_points!(
    pi::ParamInterval;
    max_recursions::Int = 2
    )

Refines profile points to make your plot more smooth. Internally uses adapted_grid to compute additional profile points. See PlotUtils.adapted_grid.

Arguments

  • pi : ParamInterval structure to update.
  • max_recursions : how many times each interval is allowed to be refined (default: 2).
source
RecipesBase.apply_recipeMethod
using Plots
plot(pi::ParamInterval)

Plots profile L(theta) for parameter theta_num, identifiability level, identifiability interval. Use update_profile_points!(pi::ProfileInterval) function to refine profile points and make your plot more smooth

source
LikelihoodProfiler.EndPointType
struct EndPoint
    value::Union{Real, Nothing}           # value of endpoint or `nothing`
    profilePoints::Array{ProfilePoint, 1} # vector of profile points
    status::Symbol                        # status
    direction::Symbol                     # `:right` or `:left`
    counter::Int                          # number of `loss_func` evaluations
    supreme::Union{Real, Nothing}         # maximal value inside profile interval
end

Structure storing one endpoint of confidence interval.

status values: :BORDER_FOUND_BY_SCAN_TOL, :BORDER_FOUND_BY_LOSS_TOL, :SCAN_BOUND_REACHED, :MAX_ITER_STOP, :LOSS_ERROR_STOP

source
LikelihoodProfiler.ParamIntervalType
struct ParamInterval
    input::ParamIntervalInput
    loss_init::Float64
    method::Symbol
    result::Tuple{EndPoint, EndPoint}
end

Structure storing result of parameter identification

source
LikelihoodProfiler.ParamIntervalInputType
struct ParamIntervalInput <: AbstractIntervalInput
    theta_init::Vector{Float64} # initial parameters vector
    theta_num::Int              # number of the parameter for identifiability analysis
    loss_func::Function         # loss function
    method::Symbol
    options::Any
end

Structure storing input data for parameter identification

source
LikelihoodProfiler.ProfilePointType
struct ProfilePoint
    value::Float64               # x value of profile point
    loss::Float64                # y value of profile point (loss function value)
    params::Array{Float64, 1}    # vector of optimal values of `loss_func` arguments
    ret::Symbol                  # `NLOpt.optimize()` return value 
    counter::Union{Int, Nothing} # number of `loss_func` evaluations
end

Structure storing one point from profile function. ret values: :FORCED_STOP, :MAXEVAL_REACHED, :FTOL_REACHED

source