Predictions
Adding a Prediction Function
In the event that a prediction function has not been added to the LikelihoodModel struct yet, we can add one using add_prediction_function!. This prediction function is for model solutions/trajectory and not the additional error we account for when predicting realisations.
LikelihoodBasedProfileWiseAnalysis.add_prediction_function! — Functionadd_prediction_function!(model::LikelihoodModel, predictfunction::Function)Adds a prediction function, predictfunction, to model and evaluates the predicted response variable(s) at the data points using the maximum likelihood estimate for model parameters. Modifies model in place.
Requirements for predictfunction
- A function to generate model predictions from that is paired with the
loglikefunction. - Takes three arguments,
θ,dataandt, in that order, whereθanddataare the same as forloglikefunctionandtneeds to be an optional third argument. - When
tis not specified, the prediction function should be evaluated for the same time points/independent variable as the data. Whentis specified, the prediction function should be evaluated for those specified time points/independent variable. A good practice is to include the time points of the data in the argument,data, so that the function can be specified aspredictfunction(θ, data, t=data.t)(heredatais aNamedTuple). - The output of the function should be a 1D vector when there is a single predicted response variable or a 2D array when there are multiple predicted response variables.
- The prediction(s) for each response variable should be stored in the columns of the array. In Julia, vectors are stored column-wise, so in the case where there is only one response variable, this will already be the case.
- The number of rows of each predicted response variable should be the same length as the vector
tused to evaluate the response.
LikelihoodBasedProfileWiseAnalysis.check_prediction_function_exists — Functioncheck_prediction_function_exists(model::LikelihoodModel)Checks if a prediction function is stored in model, returning true if there is. Otherwise false is returned and a warning is logged. Requirements for a prediction function can be seen in add_prediction_function!.
Adding an Error Function
Similarly, if a error model function (the data distribution) has not been added to the LikelihoodModel struct yet, we can add one using add_error_function!. This function is used to evaluate region reference intervals around the model trajectory given the data distribution defined as the error model.
LikelihoodBasedProfileWiseAnalysis.add_error_function! — Functionadd_error_function!(model::LikelihoodModel, errorfunction::Function)Adds an error function, errorfunction, to model. Modifies model in place. Model must contain a prediction function model.core.predictfunction to add the error function.
Requirements for errorfunction
- A function to generate lower and upper confidence quartiles (reference intervals) of each prediction realisation from
predictfunction. - Takes three arguments,
predictions,θ, andregionin that order, wherepredictionsis the array of predictions generated bypredictfunction,θis the same as forloglikefunctionandregionis the highest density region to evaluate quartiles of the error at each prediction point. - The output of the function should be a lower quartile and a upper quartile array, each with the same dimensions as
predictions.
Predefined error functions available to use
Predefined Error models
For convenience, we define example error model functions for four data distributions: Gaussian, log-normal, logit-normal and poisson. We provide versions where the standard deviation parameter σ is known (i.e. fixed at a given value) and where it's estimated (it should be a parameter within the parameter vector).
LikelihoodBasedProfileWiseAnalysis.normal_error_σ_known — Functionnormal_error_σ_known(predictions::AbstractArray,
θ::AbstractVector,
region::Float64,
σ::Real)Use a normal error model to quantify the uncertainty in the predictions of realisations, with known value of σ.
To use this function as the error model you must create a new function with only the first three arguments, which calls this function with a set value of σ.
Two equivalent examples of this specification with σ set to 1.3:
errorfunction(predictions, θ, region) = normal_error_σ_known(predictions, θ, region, 1.3)
function errorfunction(predictions, θ, region)
normal_error_σ_known(predictions, θ, region, 1.3)
endLikelihoodBasedProfileWiseAnalysis.normal_error_σ_estimated — Functionnormal_error_σ_estimated(predictions::AbstractArray,
θ::AbstractVector,
region::Float64,
σ_θindex::Int)Use a normal error model to quantify the uncertainty in the predictions of realisations, where σ is an estimated model parameter in θ.
To use this function as the error model you must create a new function with only the first three arguments, which calls this function with the index, σ_θindex, of σ in θ.
Two equivalent examples of this specification with σ stored at the end of the θ parameter vector, which has 4 elements (length 4):
errorfunction(predictions, θ, region) = normal_error_σ_estimated(predictions, θ, region, 4)
function errorfunction(predictions, θ, region)
normal_error_σ_estimated(predictions, θ, region, 4)
endLikelihoodBasedProfileWiseAnalysis.lognormal_error_σ_known — Functionlognormal_error_σ_known(predictions::AbstractArray,
θ::AbstractVector,
region::Float64,
σ::Real)Use a log normal error model to quantify the uncertainty in the predictions of realisations, with known value of σ.
To use this function as the error model you must create a new function with only the first three arguments, which calls this function with a set value of σ.
Two equivalent examples of this specification with σ set to 1.3:
errorfunction(predictions, θ, region) = lognormal_error_σ_known(predictions, θ, region, 1.3)
function errorfunction(predictions, θ, region)
lognormal_error_σ_known(predictions, θ, region, 1.3)
endLikelihoodBasedProfileWiseAnalysis.lognormal_error_σ_estimated — Functionlognormal_error_σ_estimated(predictions::AbstractArray,
θ::AbstractVector,
region::Float64,
σ_θindex::Int)Use a log normal error model to quantify the uncertainty in the predictions of realisations, where σ is an estimated model parameter in θ.
To use this function as the error model you must create a new function with only the first three arguments, which calls this function with the index, σ_θindex, of σ in θ.
Two equivalent examples of this specification with σ stored at the end of the θ parameter vector, which has 4 elements (length 4):
errorfunction(predictions, θ, region) = lognormal_error_σ_estimated(predictions, θ, region, 4)
function errorfunction(predictions, θ, region)
lognormal_error_σ_estimated(predictions, θ, region, 4)
endLikelihoodBasedProfileWiseAnalysis.logitnormal_error_σ_known — Functionlogitnormal_error_σ_known(predictions::AbstractArray,
θ::AbstractVector,
region::Float64,
σ::Real)Use a logit-normal error model to quantify the uncertainty in the predictions of realisations, with known value of σ. Predictions are required to be defined ∈ (0,1).
Output is the correct highest density region for around σ < 1.43 (at all values of predictions). If predictions are closer to 0 or 1.0 than 0.5, then higher values of σ (closer to 2) will be acceptable as the distribution will still be unimodal. Otherwise, the distribution will not be unimodal and won't identify the correct high density regions.
To use this function as the error model you must create a new function with only the first three arguments, which calls this function with a set value of σ.
Two equivalent examples of this specification with σ set to 0.9:
errorfunction(predictions, θ, region) = logitnormal_error_σ_known(predictions, θ, region, 0.9)
function errorfunction(predictions, θ, region)
logitnormal_error_σ_known(predictions, θ, region, 0.9)
endLikelihoodBasedProfileWiseAnalysis.logitnormal_error_σ_estimated — Functionlogitnormal_error_σ_estimated(predictions::AbstractArray,
θ::AbstractVector,
region::Float64,
σ_θindex::Int)Use a logit-normal error model to quantify the uncertainty in the predictions of realisations, where σ is an estimated model parameter in θ. Predictions are required to be defined ∈ (0,1).
Output is the correct highest density region for around σ < 1.43 (at all values of predictions). If predictions are closer to 0 or 1.0 than 0.5, then higher values of σ (closer to 2) will be acceptable as the distribution will still be unimodal. Otherwise, the distribution will not be unimodal and won't identify the correct high density regions.
To use this function as the error model you must create a new function with only the first three arguments, which calls this function with the index, σ_θindex, of σ in θ.
Two equivalent examples of this specification with σ stored at the end of the θ parameter vector, which has 4 elements (length 4):
errorfunction(predictions, θ, region) = logitnormal_error_σ_estimated(predictions, θ, region, 4)
function errorfunction(predictions, θ, region)
logitnormal_error_σ_estimated(predictions, θ, region, 4)
endLikelihoodBasedProfileWiseAnalysis.poisson_error — Functionpoisson_error(predictions::AbstractArray,
θ::AbstractVector,
region::Float64)Use a poisson error model to quantify the uncertainty in the predictions of realisations, where each prediction is the mean of the error model.
To use this function as the error model you don't need to create a new function, just specify this function directly.
Prediction Generation
Then to generate predictions we can use one of three functions, depending on whether we want to generate predictions from univariate or bivariate profiles, or dimensional samples.
LikelihoodBasedProfileWiseAnalysis.generate_predictions_univariate! — Functiongenerate_predictions_univariate!(model::LikelihoodModel,
t::AbstractVector,
proportion_to_keep::Real;
<keyword arguments>)Evalute and save proportion_to_keep individual predictions and their extrema from existing univariate profiles that meet the requirements of the univariate method of LikelihoodBasedProfileWiseAnalysis.desired_df_subset (see Keyword Arguments) at time points t. Modifies model in place.
Arguments
model: aLikelihoodModelcontaining model information, saved profiles and predictions.t: a vector of time points to compute predictions at.proportion_to_keep: aRealnumber ∈ [0.0,1.0] of the proportion of individual predictions to save. Default is 1.0.
Keyword Arguments
region: aRealnumber ∈ [0, 1] specifying the proportion of the density of the error model from which to evaluate the highest density region. Default is0.95.confidence_levels: a vector of confidence levels. If empty, all confidence levels of univariate profiles will be considered for evaluating predictions from. Otherwise, only confidence levels inconfidence_levelswill be considered. Default isFloat64[](any confidence level).dofs: a vector of integer degrees of freedom used to define the asymptotic threshold for the extremities of a univariate profile. If empty, all degrees of freedom for univariate profiles will be considered for evaluating predictions from. Otherwise, only degrees of freedom indofswill be considered. Default isInt[](any degree of freedom).profile_types: a vector ofAbstractProfileTypestructs. If empty, all profile types of univariate profiles are considered. Otherwise, only profiles with matching profile types will be considered. Default isAbstractProfileType[](any profile type).overwrite_predictions: boolean variable specifying whether to re-evaluate and overwrite predictions for univariate profiles that have already had predictions evaluated. Set totrueif predictions need to be evaluated for a new vector of time points. Default isfalse.show_progress: boolean variable specifying whether to display progress bars on the percentage of predictions evaluated and estimated time of completion. Default ismodel.show_progress.use_distributed: boolean variable specifying whether to use a normal for loop or a@distributedfor loop across univariate profiles. Default istrue.
Details
For each univariate profile that meets the requirements of LikelihoodBasedProfileWiseAnalysis.desired_df_subset, it uses LikelihoodBasedProfileWiseAnalysis.generate_prediction to generates the predictions for every parameter point in the profiles. The extrema of these predictions are computed (these are approximate simultaneous confidence bands for the prediction mean). The extrema and proportion_to_keep of the individual predictions are saved as a PredictionStruct in model.uni_predictions_dict, where the keys for the dictionary is the row number in model.uni_profiles_df of the corresponding profile.
Distributed Computing Implementation
If Distributed.jl is being used and use_distributed is true then the predictions from each univariate profile will be computed in parallel across Distributed.nworkers() workers.
Iteration Speed Of the Progress Meter
The time/it value is the time it takes for a prediction to be evaluated from a single point in parameter space.
LikelihoodBasedProfileWiseAnalysis.generate_predictions_bivariate! — Functiongenerate_predictions_bivariate!(model::LikelihoodModel,
t::AbstractVector,
proportion_to_keep::Real;
<keyword arguments>)Evalute and save proportion_to_keep individual predictions and their extrema from existing bivariate profiles that meet the requirements of the bivariate method of LikelihoodBasedProfileWiseAnalysis.desired_df_subset (see Keyword Arguments) at time points t. Modifies model in place.
Arguments
model: aLikelihoodModelcontaining model information, saved profiles and predictions.t: a vector of time points to compute predictions at.proportion_to_keep: aRealnumber ∈ [0.0,1.0] of the proportion of individual predictions to save. Default is 1.0.
Keyword Arguments
region: aRealnumber ∈ [0, 1] specifying the proportion of the density of the error model from which to evaluate the highest density region. Default is0.95.confidence_levels: a vector of confidence levels. If empty, all confidence levels of bivariate profiles will be considered for evaluating predictions from. Otherwise, only confidence levels inconfidence_levelswill be considered. Default isFloat64[](any confidence level).dofs: a vector of integer degrees of freedom used to define the asymptotic threshold for the boundary of a bivariate profile. If empty, all degrees of freedom for bivariate profiles will be considered for evaluating predictions from. Otherwise, only degrees of freedom indofswill be considered. Default isInt[](any degree of freedom).profile_types: a vector ofAbstractProfileTypestructs. If empty, all profile types of bivariate profiles are considered. Otherwise, only profiles with matching profile types will be considered. Default isAbstractProfileType[](any profile type).methods: a vector ofAbstractBivariateMethodstructs. If empty all methods used to find bivariate profiles are considered. Otherwise, only profiles with matching method types will be considered (struct arguments do not need to be the same). Default isAbstractBivariateMethod[](any bivariate method).overwrite_predictions: boolean variable specifying whether to re-evaluate and overwrite predictions for bivariate profiles that have already had predictions evaluated. Set totrueif predictions need to be evaluated for a new vector of time points. Default isfalse.show_progress: boolean variable specifying whether to display progress bars on the percentage of predictions evaluated and estimated time of completion. Default ismodel.show_progress.use_distributed: boolean variable specifying whether to use a normal for loop or a@distributedfor loop across bivariate profiles. Default istrue.
Details
For each bivariate profile that meets the requirements of LikelihoodBasedProfileWiseAnalysis.desired_df_subset, it uses LikelihoodBasedProfileWiseAnalysis.generate_prediction to generates the predictions for every parameter point in the profiles. The extrema of these predictions are computed (these are approximate simultaneous confidence bands for the prediction mean). The extrema and proportion_to_keep of the individual predictions are saved as a PredictionStruct in model.biv_predictions_dict, where the keys for the dictionary is the row number in model.biv_profiles_df of the corresponding profile.
Distributed Computing Implementation
If Distributed.jl is being used and use_distributed is true then the predictions from each bivariate profile will be computed in parallel across Distributed.nworkers() workers.
Iteration Speed Of the Progress Meter
The time/it value is the time it takes for a prediction to be evaluated from a single point in parameter space.
LikelihoodBasedProfileWiseAnalysis.generate_predictions_dim_samples! — Functiongenerate_predictions_dim_samples!(model::LikelihoodModel,
t::AbstractVector,
proportion_to_keep::Real;
<keyword arguments>)Evalute and save proportion_to_keep individual predictions and their extrema from existing dimensional samples that meet the requirements of the dimensional method of LikelihoodBasedProfileWiseAnalysis.desired_df_subset (see Keyword Arguments) at time points t. Modifies model in place.
Arguments
model: aLikelihoodModelcontaining model information, saved profiles and predictions.t: a vector of time points to compute predictions at.proportion_to_keep: aRealnumber ∈ [0.0,1.0] of the proportion of individual predictions to save. Default is 1.0.
Keyword Arguments
region: aRealnumber ∈ [0, 1] specifying the proportion of the density of the error model from which to evaluate the highest density region. Default is0.95.confidence_levels: a vector of confidence levels. If empty, all confidence levels of dimensional samples will be considered for evaluating predictions from. Otherwise, only confidence levels inconfidence_levelswill be considered. Default isFloat64[](any confidence level).dofs: a vector of integer degrees of freedom used to define the asymptotic threshold for the boundary of a dimensional sample. If empty, all degrees of freedom for dimensional sample will be considered for evaluating predictions from. Otherwise, only degrees of freedom indofswill be considered. Default isInt[](any degree of freedom).sample_types: a vector ofAbstractSampleTypestructs. If empty, all sample types used to find dimensional samples are considered. Otherwise, only samples with matching sample types will be considered. Default isAbstractSampleType[](any sample type).overwrite_predictions: boolean variable specifying whether to re-evaluate and overwrite predictions for dimensional samples that have already had predictions evaluated. Set totrueif predictions need to be evaluated for a new vector of time points. Default isfalse.show_progress: boolean variable specifying whether to display progress bars on the percentage of predictions evaluated and estimated time of completion. Default ismodel.show_progress.use_distributed: boolean variable specifying whether to use a normal for loop or a@distributedfor loop across dimensional samples. Default istrue.
Details
For each dimensional sample that meets the requirements of LikelihoodBasedProfileWiseAnalysis.desired_df_subset, it uses LikelihoodBasedProfileWiseAnalysis.generate_prediction to generates the predictions for every parameter point in the samples. The extrema of these predictions are computed (these are approximate simultaneous confidence bands for the prediction mean). The extrema and proportion_to_keep of the individual predictions are saved as a PredictionStruct in model.dim_predictions_dict, where the keys for the dictionary is the row number in model.dim_samples_df of the corresponding sample.
Distributed Computing Implementation
If Distributed.jl is being used and use_distributed is true then the predictions from each dimensional sample will be computed in parallel across Distributed.nworkers() workers.
Iteration Speed Of the Progress Meter
The time/it value is the time it takes for a prediction to be evaluated from a single point in parameter space.
To extract predictions give the relevant profile we can use:
LikelihoodBasedProfileWiseAnalysis.get_univariate_prediction_set — Functionget_uni_prediction_set(model::LikelihoodModel, uni_row_number::Int)Returns the PredictionStruct struct corresponding to the profile in row uni_row_number of model.uni_profiles_df.
LikelihoodBasedProfileWiseAnalysis.get_bivariate_prediction_set — Functionget_bivariate_prediction_set(model::LikelihoodModel, biv_row_number::Int)Returns the PredictionStruct struct corresponding to the profile in row biv_row_number of model.biv_profiles_df.
LikelihoodBasedProfileWiseAnalysis.get_dimensional_prediction_set — Functionget_dim_prediction_set(model::LikelihoodModel, dim_row_number::Int)Returns the PredictionStruct struct corresponding to the profile in row dim_row_number of model.dim_samples_df.
Structs
Predictions are stored in a PredictionStruct which will also contain a PredictionRealisationsStruct if the corresponding reference tolerance intervals have been evaluated.
LikelihoodBasedProfileWiseAnalysis.AbstractPredictionStruct — TypeAbstractPredictionStructSupertype for the predictions storage struct.
Subtypes
LikelihoodBasedProfileWiseAnalysis.PredictionStruct — TypePredictionStruct(predictions::Array{Real}, extrema::Array{Real}, realisations::PredictionRealisationsStruct)Struct for containing evaluated predictions corresponding to confidence profiles.
Fields
predictions: array of model predictions evaluated at the parameters given by a particular confidence profile parameter set. If a model has multiple response variables, it assumes thatmodel.core.predictfunctionstores the prediction for each variable in its columns. Values for each response variable are stored in the 3rd dimension (row=dim1, col=dim2, page/sheet=dim3). Each column corresponds to a column of the confidence profile parameter set.extrema: extrema of the predictions array.realisations: aPredictionRealisationsStructstruct.
Supertype Hiearachy
PredictionStruct <: AbstractPredictionStruct <: Any
LikelihoodBasedProfileWiseAnalysis.PredictionRealisationsStruct — TypePredictionRealisationsStruct(lq::Array{<:Real}, uq::Array{<:Real}, extrema::Array{<:Real})Struct for containing evaluated lower and upper confidence quartiles of the reference tolerance sets for prediction realisations corresponding to confidence profiles.
Fields
lq: array of the lower confidence quartile of the error model evaluated at each prediction point in a correspondingPredictionStruct.uq: array of the upper confidence quartile of the error model evaluated at each prediction point in a correspondingPredictionStruct.extrema: extrema of thelqanduqarrays.
Supertype Hiearachy
PredictionRealisationsStruct <: AbstractPredictionStruct <: Any
Index
LikelihoodBasedProfileWiseAnalysis.AbstractPredictionStructLikelihoodBasedProfileWiseAnalysis.PredictionRealisationsStructLikelihoodBasedProfileWiseAnalysis.PredictionStructLikelihoodBasedProfileWiseAnalysis.add_error_function!LikelihoodBasedProfileWiseAnalysis.add_prediction_function!LikelihoodBasedProfileWiseAnalysis.check_prediction_function_existsLikelihoodBasedProfileWiseAnalysis.generate_predictions_bivariate!LikelihoodBasedProfileWiseAnalysis.generate_predictions_dim_samples!LikelihoodBasedProfileWiseAnalysis.generate_predictions_univariate!LikelihoodBasedProfileWiseAnalysis.get_bivariate_prediction_setLikelihoodBasedProfileWiseAnalysis.get_dimensional_prediction_setLikelihoodBasedProfileWiseAnalysis.get_univariate_prediction_setLikelihoodBasedProfileWiseAnalysis.logitnormal_error_σ_estimatedLikelihoodBasedProfileWiseAnalysis.logitnormal_error_σ_knownLikelihoodBasedProfileWiseAnalysis.lognormal_error_σ_estimatedLikelihoodBasedProfileWiseAnalysis.lognormal_error_σ_knownLikelihoodBasedProfileWiseAnalysis.normal_error_σ_estimatedLikelihoodBasedProfileWiseAnalysis.normal_error_σ_knownLikelihoodBasedProfileWiseAnalysis.poisson_error