API Reference

Main entry points

ParametricAdaptiveSampling.sample_adaptive_parametricFunction
sample_adaptive_parametric(f, tmin, tmax)
sample_adaptive_parametric(f, ts_init)

Automatically sample points from a parametric function f: Float64 -> Vector of Real such that the error from linear interpolation between these points is small.

Arguments

  • f: The parameteric function. Given a Float64 in [tmin, tmax] must return an n-element tuple or vector of reals, where n must be constant across t values. This should be at least continuous on the interval.
  • tmin::Real: The lower bound of the parameter domain.
  • tmax::Real: The upper bound of the parameter domain.
  • ts_init::Collection{Real}: If this form is used, the initial t values are given explicitly. Otherwise (with tmin and tmax) we do initial sampling according to min_points. Must have size at least 2.

f must have well-defined values within [tmin, tmax], including the endpoints.

Keyword Arguments

  • errfun: Error function (Segment, Range) -> Float64 where the Range is (an estimate of) the value range. Default = err_relative_range = l-infinity distance component-wise relative to range; this is a good default for plotting.
  • tol: Tolerance to errfun for refinement. For plotting, this could be 1 / (max vertical/horizontal resolution). Default = 1e-3.
  • min_points::Integer: Minimum number of t values to sample linearly initially. Pass 2 for no presampling beyond tmin and tmax. Default = 20.
  • max_points::Real: Maximum number of points to sample. Pass Inf for no limit (though this is not recommended).

Returns

(t values :: Vector{Number}, value points :: Vector{Tuple}). For plotting applications, you'll want to discard the t values, i.e., take the second return value only.

Operation

This function works as follows: first, it performs an initial presampling step (first form) or accepts a list of initial samples (second form). Then, it considers the segment where the error of linear interpolation across the segment vs the midpoint value (midpoint in t space) is worst and splits it into two halves. Then repeat until all errors are below the tolerance or we run out of points. The estimate of the value range is updated as we add new points.

Caveats

If min_points is low, range in errfun may be catastrophically small, at least for initial refinement. You probably don't want that. More generally, if the initial samples do not represent the value range reasonably well, in some cases, refinement may add excessive points, and may even reach max_points without any meaningful progress. If this happens, you probably want to adjust initial samples.

Recursion is is by splitting into halves in t space, so if your function has details that are missed by this recursive grid, they won't show up. In this case, increasing min_points may help.

Features not yet implemented (SOMEDAY)

  • Infinite or open t ranges where we would adaptively sample larger/smaller (towards infinity or the edges) t values. This may be an instance of a more general hook into refinement.
  • Some way of detecting and avoiding polar points. Basically a way to say "it's not useful to plot this, better to cut it out". Not clear how we'd do this in a robust way without affecting some use cases negatively.
  • Optional penalty for recursion depth (or something like this) to avoid excessive concentration at polar points if max_points gets exhausted. This would have to be configured by the user with knowledge of the function, can be detrimental otherwise. Unclear if we want this.
  • Option to split into more than two parts per recursion step, or to split not into halves but by some other proportions. This could help when details are missed by the 2-split recursion (see Caveats).
  • Point density target in value space. Could help in situations where the midway point is interpolated well linearly, but there is additional variation at some other point.
  • Option to drop points that are ultimately not needed if the range expands during refinement (may reduce number of points, might be useful for some later computation steps).
  • A way to understand that the Range may actually be the wrong thing, e.g., when we use aspect_ratio=1 in the plot and the plot therefore creates additional space, or some xlim or ylim. Hard to do generically. Maybe make an option to pass the plot range explicitly.
  • Is there some smartness we can do if the (second) derivative of f is available? Choosing the splitting point better than at the interval midpoint seems attractive.
source
ParametricAdaptiveSampling.sample_adaptiveFunction
sample_adaptive(f, xmin, xmax)
sample_adaptive(f, xs)

Non-parametric variant of sample_adaptive_parametric for a function f: Float64 -> Real. See the documentation of that function. The same keyword arguments are supported.

Returns

x-y pairs :: Vector{Tuple{Real, Real}}. These are the resulting points.

Note that the return format is different from sample_adaptive_parametric because there is no separate space of "t values" here.

When this function is useful

Note: This function is often not needed. In many cases, for plotting, it will be faster to just spam a lot of x values instead. Cases where you may want this function could be:

  1. f is slow to evaluate.
  2. You want to limit the number of data points for some later computation step that is slow.
  3. You need very high precision in your output and can't affort that many points, maybe for a non-plotting use case.

In cases 2. and 3., you may also want to use custom values and/or a custom error function in the keyword arguments.

source

Other Public Members