Bivariate Profiles
The key function for evaluating bivariate profile boundaries is bivariate_confidenceprofiles!. The evaluated bivariate profile(s) will be contained within a BivariateConfidenceStruct that is stored in the LikelihoodModel.
LikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofiles! — Functionbivariate_confidenceprofiles!(model::LikelihoodModel,
θcombinations::Vector{Vector{Int}},
num_points::Int;
<keyword arguments>)Finds num_points profile_type boundary points at a specified confidence_level for each combination of two interest parameters using a specified method, optionally saving any found internal points. Saves these profiles by modifying model in place.
Arguments
model: aLikelihoodModelcontaining model information, saved profiles and predictions.θcombinations: vector of pairs of parameters to profile, as a vector of vectors of model parameter indexes.num_points: positive number of points to find on the boundary at the specified confidence level. Depending on the method, if a region of the user-provided bounds is inside the boundary some of these points will be on the bounds and inside the boundary. Set to at least 3 within the function as some methods need at least three points to work.
Keyword Arguments
confidence_level: a number ∈ (0.0, 1.0) for the confidence level on which to find theprofile_typeboundary. Default is0.95(95%).dof: an integer ∈ [2,model.core.num_pars] for the degrees of freedom used to define the asymptotic threshold (LikelihoodBasedProfileWiseAnalysis.get_target_loglikelihood) which defines the boundary of the bivariate profile. For bivariate profiles that are considered individually, it should be set to2. For profiles that are considered simultaneously, it should be set tomodel.core.num_pars. Default is2. Setting it tomodel.core.num_parsshould be reasonable when making predictions for well-identified models with<10parameters. Note: values other than2andmodel.core.num_parsmay not have a clear statistical interpretation.profile_type: whether to use the true log-likelihood function or an ellipse approximation of the log-likelihood function centred at the MLE (with optional use of parameter bounds). Available profile types areLogLikelihood,EllipseApproxandEllipseApproxAnalytical. Default isLogLikelihood()(LogLikelihood).method: a method of typeAbstractBivariateMethod. For a list of available methods usebivariate_methods()(bivariate_methods). Default isRadialRandomMethod(5)(RadialRandomMethod).θlb_nuisance: a vector of lower bounds on nuisance parameters, requireθlb_nuisance .≤ model.core.θmle. Default ismodel.core.θlb.θub_nuisance: a vector of upper bounds on nuisance parameters, requireθub_nuisance .≥ model.core.θmle. Default ismodel.core.θub.save_internal_points: boolean variable specifying whether to save points found inside the boundary during boundary computation. Internal points can be plotted in bivariate profile plots and will be used to generate predictions from a given bivariate profile. Default istrue.existing_profiles:Symbol ∈ [:ignore, :merge, :overwrite]specifying what to do if profiles already exist for a givenθcombination,confidence_level,profile_typeandmethod. See below for each symbol's meanings. Default is:merge.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. Default ismodel.find_zero_atol.optimizationsettings: aOptimizationSettingsstruct containing the optimisation settings used to find optimal values of nuisance parameters for a given pair of interest parameter values. Default ismissing(will usemodel.core.optimizationsettings).show_progress: boolean variable specifying whether to display progress bars on the percentage ofθcombinationscompleted 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 combinations of interest parameters. Set this variable tofalseif Distributed.jl is not being used. Default istrue.use_threads: boolean variable specifying, ifuse_distributedis false, whether to use parallelised for loops acrossThreads.nthreads()threads or a non-parallel for loops to find boundary points frommethodswhere boundary points are found independently. Default istrue.Fix1AxisMethodandRadialMLEMethodparallelise the finding point pair step and the finding the boundary from point pairs step.SimultaneousMethodandRadialRandomMethoddo not parallelise the finding point pair step but parallelise finding the boundary from point pairs.IterativeBoundaryMethodparallelises finding the initial boundary but not the following boundary improvement steps.AnalyticalEllipseMethoddoes not require parallelisation.
- :ignore means profiles that already exist will not be recomputed even if they contain fewer
num_pointsboundary points. - :merge means profiles that already exist will be merged with profiles from the current algorithm run to reach
num_points. If the existing profile already has at leastnum_pointsboundary points then that profile will not be recomputed. Otherwise, the specified method will be run starting from the difference betweennum_pointsand the number of points in the existing profile. The result of that method run will be merged with the existing profile. Predictions evaluated from the existing profile will be forgotten. To keep these predictions see extended help below. - :overwrite means profiles that already exist will be overwritten, regardless of how many points they contain. Predictions evaluated from the existing profile will be forgotten. To keep these predictions see extended help below.
Details
Using LikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofile this function calls the algorithm/method specified by method for each interest parameter combination in θcombinations (depending on the setting for existing_profiles and num_points if these profiles already exist). Nuisance parameters of each point in bivariate interest parameter space are found by maximising the log-likelihood function given by profile_type. Updates model.biv_profiles_df for each successful profile and saves their results as a BivariateConfidenceStruct in model.biv_profiles_dict, where the keys for the dictionary is the row number in model.biv_profiles_df of the corresponding profile. model.biv_profiles_df.num_points is the number of points found on the bivariate boundary (it does not include the number of saved internal points).
Extended help
Valid bounds
For methods that use points placed on parameter bounds to bracket for the confidence boundary, the bracketing method utilised via Roots.jl's find_zero will be unlikely to converge to the true confidence boundary for a given pair of interest parameters if the bounds on either parameter are +/- Inf or the log-likelihood function evaluates to +/- Inf. Bounds should be set to prevent this from occurring.
Preventing predictions from being forgotten when merging or overwriting profiles
To prevent predictions from being lost from existing profiles that would be overwritten when calling bivariate_confidenceprofiles!, existing profiles should be converted into a [CombinedBivariateMethod], prior to running new bivariate profiles. To do this use combine_bivariate_boundaries! on model with keyword argument not_evaluated_predictions set to false.
Distributed Computing Implementation
If Distributed.jl is being used and use_distributed is true, then the bivariate profiles of distinct interest parameter combinations will be computed in parallel across Distributed.nworkers() workers. If use_distributed is false and use_threads is true then for methods where finding boundary points is independent they will be computed in parallel across Threads.nthreads() threads for each pair of interest parameters.
Iteration Speed Of the Progress Meter
The time/it value is the time it takes for each new boundary point to be found (for all methods except for AnalyticalEllipseMethod). For AnalyticalEllipseMethod this is the time it takes to find all points on the boundary of the ellipse of two interest parameters.
bivariate_confidenceprofiles!(model::LikelihoodModel,
θcombinations_symbols::Union{Vector{Vector{Symbol}}, Vector{Tuple{Symbol, Symbol}}},
num_points::Int;
<keyword arguments>)Profiles just the provided θcombinations_symbols parameter pairs, provided as either a vector of vectors or a vector of tuples.
bivariate_confidenceprofiles!(model::LikelihoodModel,
profile_m_random_combinations::Int,
num_points::Int;
<keyword arguments>)Profiles m random two-way combinations of model parameters (sampling without replacement), where 0 < m ≤ binomial(model.core.num_pars,2).
bivariate_confidenceprofiles!(model::LikelihoodModel,
num_points::Int;
<keyword arguments>)Profiles all two-way combinations of model parameters.
LikelihoodBasedProfileWiseAnalysis.get_bivariate_confidence_set — Functionget_bivariate_confidence_set(model::LikelihoodModel, biv_row_number::Int)Returns the BivariateConfidenceStruct corresponding to the profile in row biv_row_number of model.biv_profiles_df
Methods For Finding Boundaries
We provide several heuristics for evaluating the boundaries of bivariate profiles; they're listed in order of computational efficiency. Available methods can be checked using bivariate_methods. We recommend evaluating up to 50 points; 10-30 should be appropriate for getting a reasonable approximation of a bivariate boundary, particularly if it is relatively elliptical. Note: [AnalyticalEllipseMethod] exclusively works with the [EllipseApproxAnalytical] profile type; it is also highly recommended if that profile type is of interest.
LikelihoodBasedProfileWiseAnalysis.AbstractBivariateMethod — TypeAbstractBivariateMethodSupertype for bivariate boundary finding methods. Use bivariate_methods() for a list of available methods (see bivariate_methods).
Subtypes
LikelihoodBasedProfileWiseAnalysis.AbstractBivariateVectorMethod — TypeAbstractBivariateVectorMethod <: AbstractBivariateMethodSupertype for bivariate boundary finding methods that search between two distinct points.
Subtypes
Supertype Hiearachy
AbstractBivariateVectorMethod <: AbstractBivariateMethod <: Any
LikelihoodBasedProfileWiseAnalysis.bivariate_methods — Functionbivariate_methods()Prints a list of available bivariate methods. Available bivariate methods include IterativeBoundaryMethod, RadialRandomMethod, RadialMLEMethod, SimultaneousMethod, Fix1AxisMethod and AnalyticalEllipseMethod. Note: CombinedBivariateMethod represents a destructive merge of the boundary structs of one or more methods and is not a valid method for bivariate_confidenceprofiles!.
LikelihoodBasedProfileWiseAnalysis.IterativeBoundaryMethod — TypeIterativeBoundaryMethod(initial_num_points::Int,
angle_points_per_iter::Int,
edge_points_per_iter::Int,
radial_start_point_shift::Float64=rand(),
ellipse_sqrt_distortion::Float64=1.0;
use_ellipse::Bool=false)Method for iteratively improving an initial boundary of initial_num_points, found by pushing out from the MLE point in directions defined by either RadialMLEMethod, when use_ellipse=true, or RadialRandomMethod, when use_ellipse=false (see LikelihoodBasedProfileWiseAnalysis.findNpointpairs_radialMLE!, LikelihoodBasedProfileWiseAnalysis.iterativeboundary_init and LikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofile_iterativeboundary).
Arguments
initial_num_points: a positive integer number of initial boundary points to find.angle_points_per_iter: a integer ≥ 0 for the number of edges to explore and improve based on the angle objective before switching to the edge objective. Ifangle_points_per_iter > 0andedge_points_per_iter > 0the angle objective is considered first.edge_points_per_iter: a integer ≥ 0 for the number of edges to explore and improve based on the edge objective before switching back to the angle objective. Ifangle_points_per_iter > 0andedge_points_per_iter > 0the angle objective is considered first.radial_start_point_shift: a number ∈ [0.0,1.0]. Default isrand()(defined on [0.0,1.0]), meaning that by default a different set of initial points will be found each time.ellipse_sqrt_distortion: a number ∈ [0.0,1.0]. Default is0.01.
Keyword Arguments
use_ellipse: Whether to findinitial_num_pointsby searching in directions defined by an ellipse approximation, as inRadialMLEMethod, or as inRadialRandomMethod. Default isfalse.
Details
For additional information on the radial_start_point_shift and ellipse_sqrt_distortion arguments see the keyword arguments for generate_N_clustered_points in EllipseSampling.jl.
Recommended for use with the LogLikelihood profile type. Radial directions, edge length and internal angle calculations are rescaled by the relative magnitude/scale of the two interest parameters. This is so directions and regions explored and consequently the boundary found are not dominated by the parameter with the larger magnitude.
Once an initial boundary is found by pushing out from the MLE point in directions defined by either RadialMLEMethod or RadialRandomMethod, the method seeks to improve this boundary by minimising an internal angle and an edge length objective, each considered sequentially, until the desired number of boundary points are found. As such, the method can be thought of as a mesh smoothing or improvement algorithm; we can consider the boundary found at a given moment in time to be a N-sided polygon with edges between adjacent boundary points (vertices).
Extended help
Regions we define as needing improvement are those with:
- Internal angles between adjacent edges that are far from being straight (180 degrees or π radians). The region defined by these edges is not well explored as the real boundary edge in this region likely has some degree of smooth curvature. It may also indicate that one of these edges cuts our boundary into a enclosed region and an unexplored region on the other side of the edge. In the event that a vertex is on a user-provided bound for a parameter, this objective is set to zero, as this angle region is a byproduct of user input and not the actual log-likelihood region. This objective is defined in
LikelihoodBasedProfileWiseAnalysis.internal_angle_from_pi!. - Edges between adjacent vertices that have a large euclidean distance. The regions between these vertices is not well explored as it is unknown whether the boundary actually is straight between these vertices. On average the closer two vertices are, the more likely the edge between the two points is well approximated by a straight line, and thus our mesh is a good representation of the true log-likelihood boundary. This objective is defined in
LikelihoodBasedProfileWiseAnalysis.edge_length.
The method is implemented as follows:
- Create edges between adjacent vertices on the intial boundary. Calculate angle and edge length objectives for these edges and vertices and place them into tracked heaps.
- Until found the desired number of boundary points repeat steps 3 and 4.
- For
angle_points_per_iterpoints:- Peek at the top vertex of the angle heap.
- Place a candidate point in the middle of the edge connected to this vertex that has the larger angle at the vertex the edge connects to. This is so we explore the worse defined region of the boundary.
- Use this candidate to find a new boundary point (see below).
- If found a new boundary point, break edge we placed candidate on and replace with edges to the new boundary point, updating angle and edge length objectives in the tracked heap (see
LikelihoodBasedProfileWiseAnalysis.heapupdates_success!). Else break our polygon into multiple polygons (seeLikelihoodBasedProfileWiseAnalysis.heapupdates_failure!).
- For
edge_points_per_iterpoints:- Peek at the top edge of the edge heap.
- Place a candidate point in the middle of this edge.
- Same as for step 3.3.
- Same as for step 3.4.
At least one of angle_points_per_iter and edge_points_per_iter must be non-zero.
Uses LikelihoodBasedProfileWiseAnalysis.newboundarypoint!.
If a candidate point is on the log-likelihood threshold boundary, we accept the point.
If a candidate point is inside the boundary, then we search in the normal direction to the edge until we find a boundary point or hit the parameter bound, accepting either.
If a candidate point is outside the boundary we find the edge, e_intersect of our boundary polygon that is intersected by the line in the normal direction of the candidate edge, which passes through the candidate point. Once this edge is found, we find the vertex on e_intersect that is closest to our candidate point, subject to that vertex not also being on the candidate point edge. We setup a 1D line search/bracketing method between these two points. In the event that no boundary points are found between these two points it is likely that multiple boundaries exist. If so, we break the candidate point's edge and e_intersect and reconnect the vertexes such that we now have multiple boundary polygons.
If the largest polygon has less than two points the method will display a warning message and terminate, returning the boundary found up until then.
Boundary finding method
For the initial boundary, it uses a 1D bracketing algorithm between the MLE point and the point on the user-provided bounds in the search directions to find the boundary at the desired confidence level. For boundary improvement, if the candidate point is inside the boundary, it uses a 1D bracketing algorithm between an internal point and the point on the user-provided bounds in the search direction. If the candidate point is outside the boundary, it uses a derivative-free 1D line search algorithm in the search direction. If the derivative-free approach fails, it switches to a bracketing heuristic between the candidate point and closest vertex on e_intersect.
This method is unlikely to find boundaries that do not contain the MLE point (if they exist), but it can find them. If a boundary that does not contain the MLE point is found, it is not guaranteed to be explored. In this case the the method will inform the user that multiple boundaries likely exist for this combination of model parameters.
Impact of parameter bounds
If a parameter bound is in the way of reaching the boundary in a given search direction, the point on that bound will be returned as the boundary point. For the bracketing method to work, parameter values on the bounds need to be legal and return a non Inf value from the log-likelihood function. Interest parameter bounds that have ranges magnitudes larger than the range of the boundary may prevent the true boundary from being found. Additionally, the bracketing method will be less efficient if the interest parameter bounds are far from the true boundary.
Threaded implementation
This method is partially implemented with Threads parallelisation if use_threads is set to true when calling bivariate_confidenceprofiles!. The initial boundary step is partially parallelised (the finding point pairs step is not parallelised and the finding boundary from point pairs step is parallelised) while the boundary improvement steps are not.
Internal Points
Finds between 0 and num_points - initial_num_points internal points - internal points are found when the edge being considered's midpoint is inside the boundary.
Supertype Hiearachy
IterativeBoundaryMethod <: AbstractBivariateVectorMethod <: AbstractBivariateMethod <: Any
LikelihoodBasedProfileWiseAnalysis.RadialMLEMethod — TypeRadialMLEMethod(ellipse_start_point_shift::Float64=rand(),
ellipse_sqrt_distortion::Float64=0.01)Method for finding the bivariate boundary of a confidence profile by bracketing between the MLE point and points on the provided bounds in directions given by points found on the boundary of a ellipse approximation of the log-likelihood function around the MLE, e, using EllipseSampling.jl (see LikelihoodBasedProfileWiseAnalysis.findNpointpairs_radialMLE! and LikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofile_vectorsearch).
Arguments
ellipse_start_point_shift: a number ∈ [0.0,1.0]. Default isrand()(defined on [0.0,1.0]), meaning that by default a different set of points will be found each time.ellipse_sqrt_distortion: a number ∈ [0.0,1.0]. Default is0.01.
Details
As for univariate_confidenceintervals! with keyword argument use_ellipse_approx_analytical_start=true, the profile log-likelihood function will be evaluated at each ellipse point prior to creating the point pair bracket. If:
- the ellipse point is outside of the provided bounds, the MLE point and the point on the provided bounds in the direction of the ellipse point from the MLE is the point pair for bracketing.
- the ellipse point is between the desired boundary and the provided bounds, the MLE point and the ellipse point is the point pair for bracketing.
- the ellipse point is between the MLE point and the desired boundary, the ellipse point and the point on the provided bounds is the point pair for bracketing.
- the ellipse point is exactly on the desired boundary (profile log-likelihood function evaluates to exactly the confidence level threshold), the ellipse point is the boundary point and the method will (wrongfully) state there's a point on the bounds. From a numerical perspective this is incredibly unlikely so this is regarded as not a big deal.
For additional information on arguments see the keyword arguments for generate_N_clustered_points in EllipseSampling.jl.
Recommended for use with the EllipseApprox profile type. Will produce reasonable results for the LogLikelihood profile type when bivariate profile boundaries are convex. Otherwise, IterativeBoundaryMethod, which has an option to use a starting solution from RadialMLEMethod, is preferred as it iteratively improves the quality of the boundary and can discover regions not explored by this method.
Extended help
Boundary finding method
Uses a 1D bracketing algorithm between the MLE point and the point on the user-provided bounds in the search directions to find the boundary at the desired confidence level. Will replace one of these sides with the ellipse point depending on which side of the boundary the ellipse point is on, subject to it being between the two points.
This method is unlikely to find boundaries that do not contain the MLE point (if they exist).
Impact of parameter bounds
If a parameter bound is in the way of reaching the boundary in a given search direction, the point on that bound will be returned as the boundary point. For the bracketing method to work, parameter values on the bounds need to be legal and return a non Inf value from the log-likelihood function. Interest parameter bounds that have ranges magnitudes larger than the range of the boundary may prevent the true boundary from being found. Additionally, the bracketing method will be less efficient if the interest parameter bounds are far from the true boundary.
Threaded implementation
This method is partially implemented with Threads parallelisation if use_threads is set to true when calling bivariate_confidenceprofiles!. The finding point pairs step is not parallelised and the finding boundary from point pairs step is parallelised.
Internal Points
Finds no internal points.
Supertype Hiearachy
RadialMLEMethod <: AbstractBivariateVectorMethod <: AbstractBivariateMethod <: Any
LikelihoodBasedProfileWiseAnalysis.RadialRandomMethod — TypeRadialRandomMethod(num_radial_directions::Int, use_MLE_point::Bool=true)Method for finding the bivariate boundary of a confidence profile by finding internal boundary points using a uniform random distribution on provided bounds and choosing num_radial_directions to explore from these points (see LikelihoodBasedProfileWiseAnalysis.findNpointpairs_radialrandom! and LikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofile_vectorsearch).
Arguments
num_radial_directions: an integer greater than zero.use_MLE_point: whether to use the MLE point as the first internal point found or not. Default is true (use).
Details
Recommended for use with the LogLikelihood profile type. Radial directions are rescaled by the relative magnitude/scale of the two interest parameters. This is so directions explored and consequently the boundary found are not dominated by the parameter with the larger magnitude.
The method first uniformly samples the region specified by the bounds for the two interest parameters until a point within the boundary is found. Then it chooses num_radial_directions, spaced uniformly around a circle, and rotates these directions by a random phase shift ∈ [0.0, 2π/num_directions] radians. These directions are then distorted by the relative magnitude/scale of the two interest parameters. Then for each of these directions, until it runs out of directions or finds the desired number of points, it finds the closest point on the bounds from the internal point in this direction. If the point on the bounds is inside the boundary it will be used in place of the boundary point. If the point on the bounds is outside the boundary then it will be used as the external point in the point pair. A bracketing method is then used to find a boundary point between the point pair (the bound point and the internal point). The method continues searching for internal points and generating directions until the desired number of boundary points is found.
Given a fixed number of desired boundary points, we can decrease the cost of finding internal points by increasing the number of radial directions to explore, num_radial_directions, at each internal point. However, it is important to balance num_radial_directions with the desired number of boundary points. If num_radial_directions is large relative to the number of boundary points, then the boundary the method finds may constitute a local search around the found internal points. Resultantly, there may be regions were the boundary is not well explored. This will be less of an issue for more 'circular' boundary regions.
IterativeBoundaryMethod may be preferred over this method if evaluating the log-likelihood function is expensive or the bounds provided for the interest parameters are much larger than the boundary, as the uniform random sampling for internal points will become very expensive.
Extended help
Boundary finding method
Uses a 1D bracketing algorithm between an internal point and the point on the user-provided bounds in the search direction to find the boundary at the desired confidence level.
This method can find multiple boundaries (if they exist).
Impact of parameter bounds
If a parameter bound is in the way of reaching the boundary in a given search direction, in the same way as RadialMLEMethod, the point on that bound will be returned as the boundary point. For the bracketing method to work, parameter values on the bounds need to be legal and return a non Inf value from the log-likelihood function. Interest parameter bounds that have ranges magnitudes larger than the range of the boundary may prevent the true boundary from being found. Additionally, the bracketing method will be less efficient if the interest parameter bounds are far from the true boundary.
Threaded implementation
This method is partially implemented with Threads parallelisation if use_threads is set to true when calling bivariate_confidenceprofiles!. The finding point pairs step is not parallelised and the finding boundary from point pairs step is parallelised.
Internal Points
Finds a minimum of div(num_points, num_radial_directions, RoundUp) - 1*use_MLE_point internal points.
Supertype Hiearachy
RadialRandomMethod <: AbstractBivariateVectorMethod <: AbstractBivariateMethod <: Any
LikelihoodBasedProfileWiseAnalysis.SimultaneousMethod — TypeSimultaneousMethod(min_proportion_unique::Real=0.5, include_MLE_point::Bool=true)Method for finding the bivariate boundary of a confidence profile by finding internal and external boundary points using a uniform random distribution on provided bounds, pairing these points in the order they're found and bracketing between each pair (see LikelihoodBasedProfileWiseAnalysis.findNpointpairs_simultaneous! and LikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofile_vectorsearch). Values of min_proportion_unique lower than zero may improve performance if finding either internal points or external points is difficult given the specified bounds on interest parameters.
Arguments
min_proportion_unique: a proportion ∈ (0.0, 1.0] for the minimum allowed proportion of unique points in one of the internal or external point vectors. One of these vectors will be fully unique. Whichever vector is not unique, will have at leastmin_proportion_uniqueunique points. Default is 0.5.use_MLE_point: whether to use the MLE point as the first internal point found or not. Default is true (use).
Details
Recommended for use with the LogLikelihood profile type.
When min_proportion_unique = 1.0 the method uniformly samples the region specified by the bounds for the two interest parameters until as many internal and external boundary points as the desired number of boundary points are found. These points are paired in the order they are found.
If min_proportion_unique < 1.0 the method uniformly samples the region specified by the bounds for the two interest parameters until either the internal or external number of boundary points found is as many as the number of desired boundary points. For whichever vector is less than the number of desired boundary points, we keep searching for that kind of point until at least ceil(min_proportion_unique * num_points) points have been found. Once ceil(min_proportion_unique * num_points) points have been found these points are used to fill the remainder of the vector, in order found. Internal and external points in each vector are then paired in the order they are found.
A bracketing method is then used to find a boundary point between each point pair (the external point and the internal point).
RadialRandomMethod and IterativeBoundaryMethod are preferred over this method from a computational performance standpoint as they require fewer log-likelihood evalutions (when RadialRandomMethod has parameter num_radial_directions > 1).
Extended help
Boundary finding method
Uses a 1D bracketing algorithm between a valid point pair to find the boundary at the desired confidence level.
This method can find multiple boundaries (if they exist).
Impact of parameter bounds
If a parameter bound is in the way of reaching the boundary, points will not be put on that bound. Additionally, if the true boundary is very close to a parameter bound, the method will struggle to find this region of the boundary. This is because finding the boundary in this location requires generating a random point between the boundary and the parameter bound, which becomes more difficult the closer they are.
Interest parameter bounds that have ranges magnitudes larger than the range of the boundary will make finding internal points difficult, requiring a lot of computational effort. Similarly, the inverse will be true if external points are hard to find. Smaller values of min_proportion_unique will improve performance in these cases.
The method will fail if the interest parameter bounds are fully contained by the boundary.
Threaded implementation
This method is partially implemented with Threads parallelisation if use_threads is set to true when calling bivariate_confidenceprofiles!. The finding point pairs step is not parallelised and the finding boundary from point pairs step is parallelised.
Internal Points
Finds at least ceil(min_proportion_unique*num_points) - 1*use_MLE_point internal points.
Supertype Hiearachy
SimultaneousMethod <: AbstractBivariateVectorMethod <: AbstractBivariateMethod <: Any
LikelihoodBasedProfileWiseAnalysis.Fix1AxisMethod — TypeFix1AxisMethod()Method for finding the bivariate boundary of a confidence profile by using uniform random distributions on provided bounds to draw a value for one interest parameter, fix it, and draw two values for the other interest parameter, using the pair to find a boundary point using a bracketing method if the pair are a valid bracket (see LikelihoodBasedProfileWiseAnalysis.findNpointpairs_fix1axis! and LikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofile_fix1axis).
Details
Recommended for use with the LogLikelihood profile type.
The method finds the desired number of boundary points by fixing the first and second interest parameters for half of these points each. It first draws a value for the fixed parameter using a uniform random distribution on provided bounds (e.g. ψ_x). Then it draws two values for the unfixed parameter in the same fashion (e.g. ψ_y1 and ψ_y2]). If one of these points ([ψ_x, ψ_y1] and [ψ_x, ψ_y2]) is an internal point and the other an external point, the point pair is accepted as they are a valid bracket. A bracketing method is then used to find a boundary point between the point pair (the internal and external point). The method continues searching for valid point pairs until the desired number of boundary points is found.
RadialRandomMethod and IterativeBoundaryMethod are preferred over this method from a computational performance standpoint as they require fewer log-likelihood evalutions (when RadialRandomMethod has parameter num_radial_directions > 1). SimultaneousMethod will also likely be more efficient, even though it uses four random numbers vs three, as it doesn't unneccesarily throw away points.
Extended help
Boundary finding method
Uses a 1D bracketing algorithm between a valid point pair, which are parallel to the x or y axis, to find the boundary at the desired confidence level.
This method can find multiple boundaries (if they exist).
Impact of parameter bounds
If a parameter bound is in the way of reaching the boundary, points will not be put on that bound. Additionally, if the true boundary is very close to a parameter bound, the method will struggle to find this region of the boundary. This is because finding the boundary in this location requires generating a random point between the boundary and the parameter bound, which becomes more difficult the closer they are. Interest parameter bounds that have ranges magnitudes larger than the range of the boundary will make finding internal points very difficult, requiring a lot of computational effort. Similarly, the inverse will be true if external points are hard to find. The method will fail if the interest parameter bounds are fully contained by the boundary.
Threaded implementation
This method is implemented with Threads parallelisation if use_threads is set to true when calling bivariate_confidenceprofiles!.
Internal Points
Finds num_points internal points.
Supertype Hiearachy
Fix1AxisMethod <: AbstractBivariateMethod <: Any
LikelihoodBasedProfileWiseAnalysis.AnalyticalEllipseMethod — TypeAnalyticalEllipseMethod(ellipse_start_point_shift::Float64,
ellipse_sqrt_distortion::Float64)Method for sampling the desired number of boundary points on a ellipse approximation of the log-likelihood function centred at the maximum likelihood estimate point using EllipseSampling.jl.
Arguments
ellipse_start_point_shift: a number ∈ [0.0,1.0]. Default isrand()(defined on [0.0,1.0]), meaning that by default a different set of points will be found each time.ellipse_sqrt_distortion: a number ∈ [0.0,1.0]. Default is1.0, meaning that by default points on the ellipse approximation are equally spaced with respect to arc length.
Details
Used for the EllipseApproxAnalytical profile type only: if this method is specified, then any user provided profile type will be overriden and replaced with EllipseApproxAnalytical. This ellipse approximation ignores user provided bounds.
For additional information on arguments see the keyword arguments for generate_N_clustered_points in EllipseSampling.jl.
Extended help
Boundary finding method
Explicitly finds the boundary using EllipseSampling.jl.
Internal Points
Finds no internal points.
Supertype Hiearachy
AnalyticalEllipseMethod <: AbstractBivariateMethod <: Any
Sampling Internal Points From Boundaries
In order to cheaply sample interval points within the found boundary of a bivariate profile we use sample_bivariate_internal_points!.
LikelihoodBasedProfileWiseAnalysis.sample_bivariate_internal_points! — Functionsample_bivariate_internal_points!(model::LikelihoodModel,
num_points::Int;
<keyword arguments>)Samples num_points internal points in interest parameter space of existing bivariate profiles that meet the requirements of the bivariate method of LikelihoodBasedProfileWiseAnalysis.desired_df_subset (see Keyword Arguments). Modifies model in place, with sampled internal points appended to the internal points field of each BivariateConfidenceStruct.
Arguments
model: aLikelihoodModelcontaining model information, saved profiles and predictions.num_points: number of internal points to sample from within a polygon hull approximation of a bivariate boundary and append to it's array of internal points.
Keyword Arguments
confidence_levels: a vector of confidence levels. If empty, all confidence levels of bivariate profiles will be considered for finding interval points. Otherwise, only confidence levels inconfidence_levelswill be considered. Default isFloat64[](any confidence level).dofs: a vector of integer degrees of freedom. If empty, all degrees of freedom of bivariate profiles will be considered. Otherwise, only degrees of freedom indofsare allowed. 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).sample_type: either aUniformRandomSamplesorLatinHypercubeSamplesstruct for how to sample internal points from the polygon hull.UniformRandomSamplesare homogeneously sampled from the polygon andLatinHypercubeSamplesuse the intersection of a heuristically optimised Latin Hypercube sampling plan with the polygon. Default isLatinHypercubeSamples()(LatinHypercubeSamples).hullmethod: method of typeAbstractBivariateHullMethodused to create a 2D polygon hull that approximates the bivariate boundary from a set of boundary points and internal points (method dependent). For available methods seebivariate_hull_methods(). Default isMPPHullMethod()(MPPHullMethod).θlb_nuisance: a vector of lower bounds on nuisance parameters, requireθlb_nuisance .≤ model.core.θmle. Default ismodel.core.θlb.θub_nuisance: a vector of upper bounds on nuisance parameters, requireθub_nuisance .≥ model.core.θmle. Default ismodel.core.θub.t: vector of timepoints to evaluate predictions at for each new sampled internal point from a bivariate boundary that has already had predictions evaluated. The vector must be the same vector used to produce these previous predictions, otherwise points will not be sampled from this boundary. Default ismissing.evaluate_predictions_for_samples: boolean variable specifying whether to evaluate predictions for sampled points given predictions have been evaluated for the boundary they were sampled from. Iffalse, then existing predictions will be forgotten by themodeland overwritten the next time predictions are evaluated for each profile internal points were sampled from. Default istrue.proportion_of_predictions_to_keep: The proportion of predictions fromnum_pointsinternal points to save. Does not impact the extrema calculated from predictions. Default is1.0.optimizationsettings: aOptimizationSettingsstruct containing the optimisation settings used to find optimal values of nuisance parameters for a given pair of interest parameter values. Default ismissing(will usemodel.core.optimizationsettings).show_progress: boolean variable specifying whether to display progress bars on the percentage ofθs_to_profilecompleted 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 combinations of interest parameters. Set this variable tofalseif Distributed.jl is not being used. Default istrue.use_threads: boolean variable specifying, ifuse_distributedisfalse, to use a parallelised for loop acrossThreads.nthreads()threads to evaluate the log-likelihood at each sampled point. Default istrue.
Details
For each bivariate profile that meets the requirements of LikelihoodBasedProfileWiseAnalysis.desired_df_subset it creates a 2D polygon hull from it's set of boundary and internal points (method dependent) using hullmethod and samples points from the hull using sample_type until num_points are found, rejecting any that are not inside the log-likelihood threshold at that confidence_level, dof and profile_type. For LatinHypercubeSamples this will be approximately num_points, whereas for UniformRandomSamples this will be exactly num_points. Nuisance parameters of each point in bivariate interest parameter space are found by maximising the log-likelihood function given by the profile_type of the profile.
It is highly recommended to view the docstrings of each hullmethod as the rejection rate of sampled points and the representation accuracy / coverage of the true confidence boundary varies between them, which can impact both computational performance and sampling coverage. For example, given the same set of boundary and internal points, ConvexHullMethod will produce a polygon hull that contains at least as much of the true confidence boundary as the other methods, but may have a higher rejection rate than other methods leading to higher computational cost.
Parallel Computing Implementation
If Distributed.jl is being used and use_distributed is true then the internal samples of distinct interest parameter combinations will be computed in parallel across Distributed.nworkers() workers. If use_distributed is false and use_threads is true then the internal samples of each distinct interest parameter combination will be computed in parallel across Threads.nthreads() threads. It is highly recommended to set use_threads to true in that situation.
Iteration Speed Of the Progress Meter
An iteration within the progress meter is specified as the time it takes for all internal points within a bivariate boundary to be found.
LikelihoodBasedProfileWiseAnalysis.bivariate_hull_methods — Functionbivariate_hull_methods()Prints a list of available bivariate hull methods. Available bivariate hull methods include ConvexHullMethod, ConcaveHullMethod and MPPHullMethod.
LikelihoodBasedProfileWiseAnalysis.AbstractBivariateHullMethod — TypeAbstractBivariateMethodSupertype for bivariate boundary hull methods. Use bivariate_hull_methods() for a list of available methods (see bivariate_hull_methods).
Subtypes
LikelihoodBasedProfileWiseAnalysis.ConvexHullMethod — TypeConvexHullMethod()Constructs a 2D polygon hull to sample internal points from by applying a convex hull algorithm from Meshes.jl on the collection of points given by both the boundary and saved internal points.
Details
Representation Accuracy
For convex boundaries, this method has the ability to create a more accurate representation of the true confidence boundary than MPPHullMethod, if saved internal points contain information about the boundary that is not covered by a convex hull of boundary points alone. For concave boundaries, this method will not be an accurate representation of the true confidence boundary. However, this may be desirable if the rejection rate when sampling is not too high, because the convex nature of the hull will likely contain more of the true confidence boundary than each of the other methods, particularly for lower numbers of boundary points.
Rejection Rate when Sampling Internal Points
The rejection rate when sampling internal points will be low for convex boundaries. The rejection rate may become very high for concave boundaries, if the area of the convex hull of the true confidence boundary is much larger than the area of the true confidence boundary.
Supertype Hiearachy
ConvexHullMethod <: AbstractBivariateHullMethod <: Any
LikelihoodBasedProfileWiseAnalysis.ConcaveHullMethod — TypeConcaveHullMethod()Constructs a 2D polygon hull to sample internal points from by applying a heuristic implementation of a heuristic concave hull algorithm from ConcaveHull.jl on the collection of points given by both the boundary and saved internal points (see LikelihoodBasedProfileWiseAnalysis.bivariate_concave_hull).
Details
It applies the ConcaveHull.concave_hull algorithm twice to the union of boundary and saved internal points, with the number of neighbours for the first pass chosen using a heuristic based on the number of points in the point union. Resultantly, it may result in more accurate coverage of the true confidence boundary than MPPHullMethod, for smaller numbers of boundary points, if saved internal points are in locations not enclosed by a polygon found using only boundary points. For example, if the true confidence boundary is a square and there is an internal point in the bottom left corner, but no boundary points around that corner, the boundary polygon created by this method will likely have a vertex at that internal point. However, this is not guaranteed because it is a heuristic. Bivariate methods that struggle to find boundaries close to or on the other side of a parameter bound are an example where using information on saved internal points will prove useful.
Representation Accuracy
This method has the ability to create a more accurate representation of the boundary than MPPHullMethod, but because of it's heuristic nature this is not guaranteed. However, it should be a more accurate representation of the boundary than ConvexHullMethod for non-convex boundaries (concave boundaries). The representation in the worst case where the heuristic does not return a good representation of the boundary despite a good initial point cloud, will be signficantly worse than the ConvexHullMethod and likely worse than the MPPHullMethod.
Rejection Rate when Sampling Internal Points
The rejection rate when sampling internal points will be low for convex boundaries. The rejection rate should be low for concave boundaries as well, but the nature of the heuristic used may cause concave sections to be treated as convex.
Supertype Hiearachy
ConcaveHullMethod <: AbstractBivariateHullMethod <: Any
LikelihoodBasedProfileWiseAnalysis.MPPHullMethod — TypeMPPHullMethod()Constructs a 2D polygon hull to sample internal points from by applying a minimum perimeter polygon (MPP) traveling salesman problem algorithm to the boundary points (see LikelihoodBasedProfileWiseAnalysis.minimum_perimeter_polygon!).
Details
It does not use information on the position of internal points saved while finding a boundary. This may result in less accurate coverage of the true confidence boundary for smaller numbers of boundary points. For example, if the true confidence boundary is a square and there is an internal point in the bottom left corner, but no boundary points around that corner, the boundary polygon created by this method will not enclose the area around that corner. Bivariate methods that struggle to find boundaries close to or on the other side of a parameter bound are an example where using information on saved internal points would prove useful - ConvexHullMethod or ConcaveHullMethod may be more appropriate in these cases.
Representation Accuracy
This method will create the most accurate representation of the boundary that has been found from any of the AbstractBivariateHullMethod methods, particularly for non-convex boundaries, given a sufficient number of boundary points.
Rejection Rate when Sampling Internal Points
The rejection rate when sampling internal points will be low for convex and concave boundaries.
Supertype Hiearachy
MPPHullMethod <: AbstractBivariateHullMethod <: Any
Merging Boundaries From Multiple Methods
To improve the performance of internal point sampling, it may be worth finding bivariate boundaries using a combination of methods, where one method has more guaranteed boundary coverage and the other gives a more random search of interest parameter space. For example, combining IterativeBoundaryMethod and SimultaneousMethod into a CombinedBivariateMethod.
LikelihoodBasedProfileWiseAnalysis.CombinedBivariateMethod — TypeCombinedBivariateMethod()A method representing a BivariateConfidenceStruct that has been destructively merged from one or more boundaries found with a AbstractBivariateMethod. Does not represent a method usable by bivariate_confidenceprofiles!.
LikelihoodBasedProfileWiseAnalysis.combine_bivariate_boundaries! — Functioncombine_bivariate_boundaries!(model::LikelihoodModel;
<keyword arguments>)Combines the confidence_level bivariate boundaries at dof of profile_type found using methods into a single CombinedBivariateMethod boundary for each set of interest parameters, modifying model destructively in place. Rows of model.biv_profiles_df to combine are found using the bivariate method of LikelihoodBasedProfileWiseAnalysis.desired_df_subset. Dictionary entries and dataframe rows of boundaries that have beeen combined will be deleted and the datastructures will be rebuilt according to the new row indices of model.biv_profiles_df.
Arguments
model: aLikelihoodModelcontaining model information, saved profiles and predictions.
Keyword Arguments
confidence_level: a number ∈ (0.0, 1.0) for the confidence level ofprofile_typeboundaries to combine. Default is 0.95 (95%).dof: a integer ∈ [2, model.core.numpars] for the degrees of freedom of `profiletype` boundaries to combine. Default is 2.profile_type: the profile type of boundaries to combine. Default isLogLikelihood()(LogLikelihood).methods: a vector of methods of typeAbstractBivariateMethodfor combining boundaries found using those method types.methodsshould not containCombinedBivariateMethod, but the case where it is included inmethodsis handled: it will be removed from the vector. Default isAbstractBivariateMethod[](boundaries found using all methods are combined).not_evaluated_predictions: a boolean specifiying whether to combine only boundaries that have not had or have had predictions evaluated. If predictions are evaluated for the combined struct (if it exists) but not for the rows to combine with it, they will not be combined, and vice versa. Default istrue(combine boundaries that have not had predictions evaluated).
If predictions have been evaluated: the time points at which predictions have been evaluated at must be the same for all of the boundaires that are being combined.
Index
LikelihoodBasedProfileWiseAnalysis.AbstractBivariateHullMethodLikelihoodBasedProfileWiseAnalysis.AbstractBivariateMethodLikelihoodBasedProfileWiseAnalysis.AbstractBivariateVectorMethodLikelihoodBasedProfileWiseAnalysis.AnalyticalEllipseMethodLikelihoodBasedProfileWiseAnalysis.CombinedBivariateMethodLikelihoodBasedProfileWiseAnalysis.ConcaveHullMethodLikelihoodBasedProfileWiseAnalysis.ConvexHullMethodLikelihoodBasedProfileWiseAnalysis.Fix1AxisMethodLikelihoodBasedProfileWiseAnalysis.IterativeBoundaryMethodLikelihoodBasedProfileWiseAnalysis.MPPHullMethodLikelihoodBasedProfileWiseAnalysis.RadialMLEMethodLikelihoodBasedProfileWiseAnalysis.RadialRandomMethodLikelihoodBasedProfileWiseAnalysis.SimultaneousMethodLikelihoodBasedProfileWiseAnalysis.bivariate_confidenceprofiles!LikelihoodBasedProfileWiseAnalysis.bivariate_hull_methodsLikelihoodBasedProfileWiseAnalysis.bivariate_methodsLikelihoodBasedProfileWiseAnalysis.combine_bivariate_boundaries!LikelihoodBasedProfileWiseAnalysis.get_bivariate_confidence_setLikelihoodBasedProfileWiseAnalysis.sample_bivariate_internal_points!