Initialisation
Model Initialisation
To initialise a model for profile likelihood evaluation we use initialise_LikelihoodModel which returns a struct of type LikelihoodModel. This struct contains all the information we require for the PWA workflow and will also contain computed profiles and profile-wise predictions.
LikelihoodBasedProfileWiseAnalysis.initialise_LikelihoodModel — Functioninitialise_LikelihoodModel(loglikefunction::Function,
predictfunction::Union{Function, Missing},
errorfunction::Union{Function, Missing}
data::Union{Tuple, NamedTuple},
θnams::Vector{<:Symbol},
θinitialguess::AbstractVector{<:Real},
θlb::AbstractVector{<:Real},
θub::AbstractVector{<:Real},
θmagnitudes::AbstractVector{<:Real}=Float64[];
<keyword arguments>)Initialises a LikelihoodModel struct, which contains all model information, profiles, samples and predictions. Solves for the maximum likelihood estimate of loglikefunction.
Arguments
loglikefunction: a log-likelihood function to maximise which takes two arguments,θanddata, in that order, where θ is a vector containing the values of each parameter inθnamesanddatais a Tuple or NamedTuple - seedatabelow. Set up to be used in a maximisation objective.predictfunction: a prediction function to generate model predictions from that is paired with theloglikefunction. Requirements for the prediction function can be seen inadd_prediction_function!. It can also bemissingif no function is provided toinitialise_LikelihoodModel, because predictions are not required when evaluating parameter profiles. The function can be added at a later point usingadd_prediction_function!.errorfunction: an error function used to predict realisations from predictions generated withpredictfunction. Requirements for the error function can be seen inadd_error_function!. It can also bemissingif no function is provided toinitialise_LikelihoodModel, because predictions are not required when evaluating parameter profiles. The function can be added at a later point usingadd_error_function!.data: a Tuple or a NamedTuple containing any additional information required by the log-likelihood function, such as the time points to be evaluated at.θnames: a vector of symbols containing the names of each parameter, e.g.[:λ, :K, :C0].θinitialguess: a vector containing the initial guess for the values of each parameter. Used to find the MLE point.θlb: a vector of lower bounds on parameters.θub: a vector of upper bounds on parameters.θmagnitudes: a vector of the relative magnitude of each parameter. If not provided, it will be estimated using the difference ofθlbandθubwithLikelihoodBasedProfileWiseAnalysis.calculate_θmagnitudes. Can be updated after initialisation usingsetmagnitudes!.
Keyword Arguments
optimizationsettings: optimization settings used to optimize the log-likelihood function using Optimization.jl. Default isdefault_OptimizationSettings()(seedefault_OptimizationSettings).uni_row_prealloaction_size: number of rows ofuni_profiles_dfto preallocate. Default ismissing(a single row).biv_row_preallocation_size: number of rows ofbiv_profiles_dfto preallocate. Default ismissing(a single row).dim_row_preallocation_size: number of rows ofdim_samples_dfto preallocate. Default ismissing(a single row).find_zero_atol: aRealnumber greater than zero for the absolute tolerance of the log-likelihood function value from the target value to be used when searching for confidence intervals/boundaries. Default is0.001.show_progress: Whether to show the progress of profiling and predictions.
If you initialise an array within the provided log-likelihood function, e.g. using zeros, then for automatic differentiation methods to work you need to also initialise the type of the array to be based on the type of the input values of θ. Otherwise, zeros will by default create an array with element types Float64 which will likely return errors. For example, use:
my_new_array = zeros(eltype(θ), dimensions)where eltype passes the element type of the θ vector into the zeros function.
initialise_LikelihoodModel(loglikefunction::Function,
predictfunction::Function,
data::Union{Tuple, NamedTuple},
θnames::Vector{<:Symbol},
θinitialGuess::Vector{<:Float64},
θlb::Vector{<:Float64},
θub::Vector{<:Float64},
θmagnitudes::Vector{<:Real}=zeros(0);
<keyword arguments>)Alternate version of initialise_LikelihoodModel that can be called without a error function. The function can be added at a later point using add_error_function!.
initialise_LikelihoodModel(loglikefunction::Function,
data::Union{Tuple, NamedTuple},
θnames::Vector{<:Symbol},
θinitialGuess::Vector{<:Float64},
θlb::Vector{<:Float64},
θub::Vector{<:Float64},
θmagnitudes::Vector{<:Real}=zeros(0);
<keyword arguments>)Alternate version of initialise_LikelihoodModel that can be called without a prediction and error function. The functions can be added at a later point using add_prediction_function! and add_error_function!.
Model Representation
LikelihoodBasedProfileWiseAnalysis.LikelihoodModel — TypeLikelihoodModel(core::Union{CoreLikelihoodModel, BaseLikelihoodModel},
ellipse_MLE_approx::Union{Missing, EllipseMLEApprox},
num_uni_profiles::Int,
num_biv_profiles::Int,
num_dim_samples::Int,
uni_profiles_df::DataFrame,
biv_profiles_df::DataFrame,
dim_samples_df::DataFrame,
uni_profile_row_exists::Dict{Tuple{Int, AbstractProfileType}, DefaultDict{Float64, Int}},
biv_profile_row_exists::Dict{Tuple{Tuple{Int, Int}, AbstractProfileType, AbstractBivariateMethod}, DefaultDict{Float64, Int}},
dim_samples_row_exists::Dict{Union{AbstractSampleType, Tuple{Vector{Int}, AbstractSampleType}}, DefaultDict{Float64, Int}},
uni_profiles_dict::Dict{Int, UnivariateConfidenceStruct},
biv_profiles_dict::Dict{Int, BivariateConfidenceStruct},
dim_samples_dict::Dict{Int, SampledConfidenceStruct},
uni_predictions_dict::Dict{Int, AbstractPredictionStruct},
biv_predictions_dict::Dict{Int, AbstractPredictionStruct},
dim_predictions_dict::Dict{Int, AbstractPredictionStruct},
show_progress::Bool)Mutable struct containing all the information required to compute profiles, samples and predictions. Created by initialise_LikelihoodModel.
Fields
core: aCoreLikelihoodModelorBaseLikelihoodModelstruct.ellipse_MLE_approx: aEllipseMLEApproxstruct OR a missing value if the ellipse approximation of the log-likelihood at the MLE point has not been evaluated yet.num_uni_profiles: the number of different univariate profiles that have been evaluated (distinct combinations of different confidence levels,AbstractProfileTypestructs and single interest parameters). Specifies the number of valid rows inuni_profiles_df.num_biv_profiles: the number of different bivariate profiles that have been evaluated (distinct combinations of different confidence levels,AbstractProfileTypestructs,AbstractBivariateMethodstructs and two interest parameters). Specifies the number of valid rows inbiv_profiles_df.num_dim_samples: the number of different dimensional profiles that have been evaluated (distinct combinations of different confidence levels,AbstractProfileTypestructs,AbstractSampleTypestructs and sets of interest parameters). Specifies the number of valid rows indim_samples_df.uni_profiles_df: a DataFrame with each row containing information on each univariate profile evaluated, where the row index is the key for that profile inuni_profiles_dictanduni_predictions_dict.biv_profiles_df: a DataFrame with each row containing information on each bivariate profile evaluated, where the row index is the key for that profile inbiv_profiles_dictandbiv_predictions_dict.dim_samples_df: a DataFrame with each row containing information on each dimensional sample evaluated, where the row index is the key for that sample indim_samples_dictanddim_predictions_dict.uni_profile_row_exists: a dictionary containing information on whether a row inuni_profiles_dfexists for a given combination of interest parameter, degrees of freedom,AbstractProfileTypeand confidence level. If it does exist, the value ofuni_profile_row_exists[(θi, dof, profile_type)][confidence_level]will be the row index inuni_profiles_df, otherwise it will be0.biv_profile_row_exists: a dictionary containing information on whether a row inbiv_profiles_dfexists for a given combination of two interest parameters,AbstractProfileType,AbstractBivariateMethodand confidence level. If it does exist, the value ofbiv_profile_row_exists[((θi, θj), dof, profile_type, bivariate_method)][confidence_level]will be the row index inbiv_profiles_dfotherwise it will be0.dim_samples_row_exists: a dictionary containing information on whether a row indim_samples_dfexists for a given combination of interest parameter,AbstractProfileType,AbstractSampleTypeand confidence level. If it does exist, it's value will be the row index indim_samples_dfotherwise it will be0. For a full likelihood sample this can be accessed usingdim_samples_row_exists[sample_type][confidence_level]. For a lower dimensional sample this can be accessed usingdim_samples_row_exists[(θindices, dof, sample_type)][confidence_level].uni_profiles_dict: a dictionary with keys of type Integer and values of typeUnivariateConfidenceStructcontaining the profile for each valid row inuni_profiles_df. The row index ofuni_profiles_dfis the key for the corresponding profile.biv_profiles_dict: a dictionary with keys of type Integer and values of typeBivariateConfidenceStructcontaining the profile for each valid row inbiv_profiles_df. The row index ofbiv_profiles_dfis the key for the corresponding profile.dim_samples_dict: a dictionary with keys of type Integer and values of typeSampledConfidenceStructcontaining the profile for each valid row indim_samples_df. The row index ofdim_samples_dfis the key for the corresponding profile.uni_predictions_dict: a dictionary with keys of type Integer and values of typePredictionStructcontaining the predictions from the profiles inuni_profiles_dictfor each valid row inuni_profiles_df. The row index ofuni_profiles_dfis the key for the corresponding prediction, if that prediction has been calculated usinggenerate_predictions_univariate!.biv_predictions_dict: a dictionary with keys of type Integer and values of typePredictionStructcontaining the predictions from the profiles inbiv_profiles_dictfor each valid row inbiv_profiles_df. The row index ofbiv_profiles_dfis the key for the corresponding prediction, if that prediction has been calculated usinggenerate_predictions_bivariate!.dim_predictions_dict: a dictionary with keys of type Integer and values of typePredictionStructcontaining the predictions from the profiles indim_samples_dictfor each valid row indim_samples_df. The row index ofdim_samples_dfis the key for the corresponding prediction, if that prediction has been calculated usinggenerate_predictions_dim_samples!.find_zero_atol: aRealnumber greater than zero for the absolute tolerance of the log-likelihood function value from the target value to be used when searching for confidence intervals/boundaries.show_progress: a boolean specifying whether to show the progress of profile methods with respect to sets of interest parameter(s).
Supertype Hiearachy
LikelihoodModel <: Any
LikelihoodBasedProfileWiseAnalysis.CoreLikelihoodModel — TypeCoreLikelihoodModel(loglikefunction::Function,
predictfunction::Union{Function, Missing},
errorfunction::Union{Function, Missing}
data::Union{Tuple, NamedTuple},
θnames::AbstractVector
θname_to_index::Dict{Symbol, Int},
θlb::AbstractVector{<:Real},
θub::AbstractVector{<:Real},
θmagnitudes::AbstractVector{<:Real},
θmle::Vector{<:Float64},
ymle::Array{<:Real},
maximisedmle::Float64,
num_pars::Int)Struct containing the core information required to define a LikelihoodModel. For additional information on parameters (where repeated), see initialise_LikelihoodModel.
Fields
loglikefunction: a log-likelihood function which takes two arguments,θanddata, in that order. Set up to be used in a maximisation objective.predictfunction: a prediction function to generate model predictions from that is paired with theloglikefunction.errorfunction: an error function used to predict realisations from predictions generated withpredictfunction.optimizationsettings: aOptimizationSettingsstruct of settings to use for optimisation unless others are specified.data: a Tuple or a NamedTuple containing any additional information required by the log-likelihood function, such as the time points to be evaluated at.θnames: a vector of symbols containing the names of each parameter, e.g.[:λ, :K, :C0].θname_to_index: a dictionary with keys of type Symbol and values of type Int, with the key being an element ofθnamesand the value being the corresponding index of the key inθnames.θlb: a vector of lower bounds on parameters.θub: a vector of upper bounds on parameters.θmagnitudes: a vector of the relative magnitude of each parameter.θmle: a vector containing the maximum likelihood estimate for each parameter.ymle: an array containing the output of the prediction function atθmleanddata.maximisedmle: the value of the log-likelihood function evaluated atθmle.num_pars: the number of model parameters,length(θnames).
Supertype Hiearachy
CoreLikelihoodModel <: Any
LikelihoodBasedProfileWiseAnalysis.BaseLikelihoodModel — TypeBaseLikelihoodModel(data::Union{Tuple, NamedTuple},
θnames::AbstractVector,
θname_to_index::Dict{Symbol, Int},
θlb::AbstractVector{<:Real},
θub::AbstractVector{<:Real},
θmagnitudes::AbstractVector{<:Real},
θmle::Vector{<:Float64},
ymle::Array{<:Real},
maximisedmle::Float64,
num_pars::Int)Version of CoreLikelihoodModel without functions, allowing loading of the model struct without those functions defined.
Supertype Hiearachy
BaseLikelihoodModel <: Any
Ellipse Approximation
In order to evaluate the ellipse approximation of the normalised log-likelihood function we need to evaluate the observed Fisher information matrix (FIM) [2]. The observed FIM is a quadratic approximation of the curvature of the log-likelihood function at the maximum likelihood estimate (MLE) for parameters. This then allows us to define the ellipse approximation of the log-likelihood function, represented as the EllipseApproxAnalytical and EllipseApprox profile types. These have corresponding equations to evaluate their approximation to the log-likelihood function, as given by LikelihoodBasedProfileWiseAnalysis.analytic_ellipse_loglike and LikelihoodBasedProfileWiseAnalysis.ellipse_loglike, respectively.
LikelihoodBasedProfileWiseAnalysis.getMLE_ellipse_approximation! — FunctiongetMLE_ellipse_approximation!(model::LikelihoodModel)Creates the ellipse approximation of the log-likelihood function at the maximum likelihood estimate, modifying model in place, computing the negative hessian of the log-likelihood function and it's pseudoinverse using getMLE_hessian_and_covariance. These matrices are stored as a EllipseMLEApprox struct within model at model.ellipse_MLE_approx.
LikelihoodBasedProfileWiseAnalysis.check_ellipse_approx_exists! — Functioncheck_ellipse_approx_exists!(model::LikelihoodModel)Checks if the ellipse approximation at the maximum likelihood estimate has been created and if not creates it using getMLE_ellipse_approximation!, modifying model in place.
LikelihoodBasedProfileWiseAnalysis.EllipseMLEApprox — TypeEllipseMLEApprox(Hmle::Matrix{<:Float64}, Γmle::Matrix{<:Float64})Struct containing two n*n arrays representing the ellipse approximation of the log-likelihood function around the MLE point. See getMLE_ellipse_approximation!
Fields
Hmle: a n*n array, where n is the number of model parameters, containing the negative Hessian of the log-likelihood function evaluated at the MLE point.Γmle: a n*n array, where n is the number of model parameters, containing the inverse ofHmle.
Supertype Hiearachy
EllipseMLEApprox <: Any
Modifying Parameter Magnitudes and Bounds
If desired, the defined parameter magnitudes and bounds contained with a LikelihoodModel can be updated using setmagnitudes! and setbounds!.
LikelihoodBasedProfileWiseAnalysis.setmagnitudes! — Functionsetmagnitudes!(model::LikelihoodModel, θmagnitudes::AbstractVector{<:Real})Updates the magnitudes of each parameter in model from model.core.θmagnitudes to θmagnitudes.
LikelihoodBasedProfileWiseAnalysis.setbounds! — Functionsetbounds!(model::LikelihoodModel;
lb::AbstractVector{<:Real}=Float64[],
ub::AbstractVector{<:Real}=Float64[])Updates the parameter bounds in model from model.core.θlb to lb if specified and from model.core.θub to ub if specified. lb and ub are keyword arguments.
Parameter Transformations
To assist with parameter transformations we provide functions for transforming parameter bounds given a monotonic forward mapping from the original parameter space to the transformed parameter space.
LikelihoodBasedProfileWiseAnalysis.transformbounds — Functiontransformbounds(transformfun::Function,
lb::AbstractVector{<:Real},
ub::AbstractVector{<:Real},
independentParameterIndexes::Vector{<:Int}=Int[],
dependentParameterIndexes::Vector{<:Int}=Int[])Given a monotonic (increasing or decreasing) function, transformfun, that describes a parameter transformation, return the lower and upper bounds in the transformed space that correspond to the lower and upper bounds in the original space. Uses a heuristic to evaluate the bound transformation. We assume that the ordering of parameters stay the same for the purposes of independentParameterIndexes and dependentParameterIndexes.
Arguments
transformfun: a function describing the forward transformation between parameter spaceθandΘ. Should take in a single argument,θ, a vector of parameter values in the original space and returnΘ, a vector parameter values in the transformed space. These vectors need to be the same length aslbandub.lb: a vector of lower bounds on parameters.ub: a vector of upper bounds on parameters.independentParameterIndexes: a vector of parameter indexes where the new parameterΘ[i]depends only ontransformfun(θ[i]).dependentParameterIndexes: a vector of parameter indexes where the new parameterΘ[i]depends ontransformfun(θ[i], θ[j], j!=i).
The heuristic used for dependent parameters may fail if there are multiple local minima for the appropriate bounds to use. In this case transformbounds_NLopt should be used.
Warns if any of the returned bounds are Inf using LikelihoodBasedProfileWiseAnalysis.checkforInf.
LikelihoodBasedProfileWiseAnalysis.transformbounds_NLopt — Functiontransformbounds_NLopt(transformfun::Function,
lb::AbstractVector{<:Real},
ub::AbstractVector{<:Real})Given a monotonic (increasing or decreasing) function, transformfun, that describes a parameter transformation, return the lower and upper bounds in the transformed space that correspond to the lower and upper bounds in the original space. Uses a naturally binary integer programme if transformfun is monotonic on θ between lb and ub.
Arguments
transformfun: a function describing the forward transformation between parameter spaceθandΘ. Should take in a single argument,θ, a vector of parameter values in the original space and returnΘ, a vector parameter values in the transformed space. These vectors need to be the same length aslbandub.lb: a vector of lower bounds on parameters.ub: a vector of upper bounds on parameters.
Warns if any of the returned bounds are Inf using LikelihoodBasedProfileWiseAnalysis.checkforInf.
Example Usage
For example, if we want to log transform the parameterisation given a defined log-likelihood function and lower and upper bounds we:
Calculate the new, log-transformed bounds using the forward transformation log.(θ). This is easy enough to be done without transformbounds_NLopt, but other more complex transformations may not be.
lb, ub = ..., ...
f(x) = log.(x)
lb_log, ub_log = transformbounds_NLopt(f, lb, ub)Define a new log-likelihood function which uses the backward transformation exp.(θ) from the transformed parameter space to the original space.
function loglike(θ, data)
...
end
function loglike_log(Θ, data)
return loglike(exp.(Θ), data)
endOptimization Settings
We can set our default optimisation settings using a OptimizationSettings struct. This will be contained within the CoreLikelihoodModel field of a LikelihoodModel and can be passed as an option to initialise_LikelihoodModel. Unless different ones are passed to functions for computing profiles, they will also be used for that purpose. It may be useful to compute the MLE using conservative settings for accuracy, and then use less conservative settings for the optimisation of nuisance parameters along parameter profiles.
LikelihoodBasedProfileWiseAnalysis.default_OptimizationSettings — Functiondefault_OptimizationSettings()Creates a OptimizationSettings struct with defaults of:
adtype:SciMLBase.NoAD()(no automatic differentiation).solve_alg:NLopt.LN_BOBYQA().solve_kwargs:(maxtime=15, xtol_rel=1e-12).
If this function causes an error then LikelihoodBasedProfileWiseAnalysis needs to be loaded. Alternatively, the packages SciMLBase, Optimization and OptimizationNLopt need to be loaded.
LikelihoodBasedProfileWiseAnalysis.create_OptimizationSettings — Functioncreate_OptimizationSettings(model::LikelihoodModel;
adtype::Union{SciMLBase.AbstractADType, Missing}=missing,
solve_alg=missing,
solve_kwargs::Union{NamedTuple, Missing}=missing)Method for creating a OptimizationSettings struct with each field of the struct as a keyword argument. If a keyword argument is not provided, then the setting in model.core.optimizationsettings is used (the currently set optimization settings).
create_OptimizationSettings(;
adtype::Union{SciMLBase.AbstractADType, Missing}=missing,
solve_alg=missing,
solve_kwargs::Union{NamedTuple, Missing}=missing)Method for creating a OptimizationSettings struct with each field of the struct as a keyword argument. If a keyword argument is not provided, then the default setting in default_OptimizationSettings is used.
LikelihoodBasedProfileWiseAnalysis.set_OptimizationSettings! — Functionset_OptimizationSettings!(model::LikelihoodModel, optimizationsettings::OptimizationSettings)Updates the optimization settings contained with model with optimizationsettings.
LikelihoodBasedProfileWiseAnalysis.OptimizationSettings — TypeOptimizationSettings(adtype::SciMLBase.AbstractADType,
solve_alg,
solve_kwargs::NamedTuple)Contains the optimization settings used to solve for nuisance parameters of the log-likelihood function using Optimization.jl. Note: Optimization.jl uses a minimisation objective, while the log-likelihood functions are setup with a maximisation objective in mind. This is handled inside the package. Defaults can be seen in default_OptimizationSettings.
Fields
adtype: a method for automatically generating derivative functions of the log-likelihood function being optimised. For available settings see Automatic Differentiation Construction Choice Recommendations. Note: the corresponding package toadtypeneeds to be loaded withusing. e.g. settingadtype = AutoFiniteDiff()requiresusing FiniteDiff. Derivative-based algorithms insolve_algwill require anadtypeto be specified.AutoFiniteDiff(),AutoForwardDiff(),AutoReverseDiff(),AutoZygote()andAutoTracker()have been tested to work for finding the maximum likelihood estimate with regular model parameters.- For profiling,
AutoFiniteDiff(),AutoForwardDiff(),AutoReverseDiff()andAutoTracker()have been tested to work with regular model parameters. - If the variance of the error distribution is included as a parameter to be estimated,
AutoFiniteDiff(),AutoForwardDiff()andAutoTracker()have been tested to work for finding the MLE and profiling. AutoFiniteDiff()will always work, regardless of model specification, but may be less optimal than other methods.
solve_alg: an algorithm to use to solve for the nuisance parameters of the log-likelihood function defined within Optimization.jl. The package is loaded with the Optimization integration of NLopt, so any of the NLopt algorithms are available without having to load another package (see OptimizationNLopt). Good starting methods may beNLopt.LN_BOBYQA(),NLopt.LN_NELDERMEAD()andNLopt.LD_LBFGS. Other packages can be used as well - see Overview of the Optimizers.solve_kwargs: aNamedTupleof keyword arguments used to set solver options likemaxitersandmaxtime. For a list of common solver arguments see: Common Solver Options. Other specific package arguments may also be available. For NLopt see Methods.
Supertype Hiearachy
OptimizationSettings <: Any
LikelihoodBasedProfileWiseAnalysis.optimise — Functionoptimise(fun, q, options::OptimizationSettings, θ₀, lb, ub)Optimization.jl optimiser used for calculating the values of nuisance parameters. Default values of options use NLopt.jl algorithms (seedefault_OptimizationSettings).
optimise(fun, options::OptimizationSettings, θ₀, lb, ub)Optimization.jl optimiser used for calculating the bound transformations. Default values of options use NLopt.jl algorithms (seedefault_OptimizationSettings).
LikelihoodBasedProfileWiseAnalysis.optimise_unbounded — Functionoptimise_unbounded(fun, q, options::OptimizationSettings, θ₀)Alternative version of optimise without nuisance parameter bounds. Used for computing the nuisance parameters of EllipseApproxAnalytical profiles. Default values of options use NLopt.jl algorithms (seedefault_OptimizationSettings).
Index
LikelihoodBasedProfileWiseAnalysis.BaseLikelihoodModelLikelihoodBasedProfileWiseAnalysis.CoreLikelihoodModelLikelihoodBasedProfileWiseAnalysis.EllipseMLEApproxLikelihoodBasedProfileWiseAnalysis.LikelihoodModelLikelihoodBasedProfileWiseAnalysis.OptimizationSettingsLikelihoodBasedProfileWiseAnalysis.check_ellipse_approx_exists!LikelihoodBasedProfileWiseAnalysis.create_OptimizationSettingsLikelihoodBasedProfileWiseAnalysis.default_OptimizationSettingsLikelihoodBasedProfileWiseAnalysis.getMLE_ellipse_approximation!LikelihoodBasedProfileWiseAnalysis.initialise_LikelihoodModelLikelihoodBasedProfileWiseAnalysis.optimiseLikelihoodBasedProfileWiseAnalysis.optimise_unboundedLikelihoodBasedProfileWiseAnalysis.set_OptimizationSettings!LikelihoodBasedProfileWiseAnalysis.setbounds!LikelihoodBasedProfileWiseAnalysis.setmagnitudes!LikelihoodBasedProfileWiseAnalysis.transformboundsLikelihoodBasedProfileWiseAnalysis.transformbounds_NLopt