Gaussian Processes Reference

GpABC functions for Gaussian Process Regression. See also Gaussian Processes Overview, Gaussian Processes Examples.

Index

Types and Functions

GpABC.GPModelType
GPModel(training_x::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}},
        training_y::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}},
        kernel::AbstractGPKernel
        [,test_x::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}}=zeros(0,0)])

Constructor of GPModel that allows the kernel to be specified. Arguments that are passed as 1-d vectors will be reshaped into 2-d.

source
GpABC.GPModelType
GPModel(training_x::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}},
        training_y::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}}
        [,test_x::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}}=zeros(0,0)])

Default constructor of GPModel, that will use SquaredExponentialIsoKernel. Arguments that are passed as 1-d vectors will be reshaped into 2-d.

source
GpABC.GPModelType
GPModel

The main type that is used by most functions within the package.

All data matrices are row-major.

Fields

  • kernel::AbstractGPKernel: the kernel
  • gp_training_x::AbstractArray{Float64, 2}: training x. Size: $n \times d$.
  • gp_training_y::AbstractArray{Float64, 2}: training y. Size: $n \times 1$.
  • gp_test_x::AbstractArray{Float64, 2}: test x. Size: $m \times d$.
  • gp_hyperparameters::AbstractArray{Float64, 1}: kernel hyperparameters, followed by standard deviation of intrinsic noise $\sigma_n$, which is always the last element in the array.
  • cache::HPOptimisationCache: cache of matrices that can be re-used between calls to gp_loglikelihood and gp_loglikelihood_grad
source
GpABC.GPModelMethod
GPModel(;training_x::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}}=zeros(0,0),
    training_y::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}}=zeros(0,0),
    test_x::Union{AbstractArray{Float64, 2}, AbstractArray{Float64, 1}}=zeros(0,0),
    kernel::AbstractGPKernel=SquaredExponentialIsoKernel(),
    gp_hyperparameters::AbstractArray{Float64, 1}=Array{Float64}(0))

Constructor of GPModel with explicit arguments. Arguments that are passed as 1-d vectors will be reshaped into 2-d.

source
GpABC.gp_loglikelihoodMethod
gp_loglikelihood(gpm::GPModel)

Compute the log likelihood function, based on the kernel and training data specified in gpm.

\[log p(y \vert X, \theta) = - \frac{1}{2}(y^TK^{-1}y + log \vert K \vert + n log 2 \pi)\]

source
GpABC.gp_loglikelihood_logMethod
gp_loglikelihood_log(theta::AbstractArray{Float64, 1}, gpm::GPModel)

Log likelihood function with log hyperparameters. This is the target function of the hyperparameters optimisation procedure. Its gradient is computed by gp_loglikelihood_grad.

source
GpABC.gp_regressionMethod
gp_regression(test_x::Union{AbstractArray{Float64, 1}, AbstractArray{Float64, 2}},
    gpem::GPModel; <optional keyword arguments>)
source
GpABC.gp_regressionMethod
gp_regression(gpm::GPModel; <optional keyword arguments>)

Run the Gaussian Process Regression.

Arguments

  • gpm: the GPModel, that contains the training data (x and y), the kernel, the hyperparameters and the test data for running the regression.
  • test_x: if specified, overrides the test x in gpm. Size $m \times d$.
  • log_level::Int (optional): log level. Default is 0, which is no logging at all. 1 makes gp_regression print basic information to standard output.
  • full_covariance_matrix::Bool (optional): whether we need the full covariance matrix, or just the variance vector. Defaults to false (i.e. just the variance).
  • batch_size::Int (optional): If full_covariance_matrix is set to false, then the mean and variance vectors will be computed in batches of this size, to avoid allocating huge matrices. Defaults to 1000.
  • observation_noise::Bool (optional): whether the observation noise (with variance $\sigma_n^2$) should be included in the output variance. Defaults to true.

Return

A tuple of (mean, var). mean is a mean vector of the output multivariate Normal distribution, of size $m$. var is either the covariance matrix of size $m \times m$, or a variance vector of size $m$, depending on full_covariance_matrix flag.

source
GpABC.gp_regression_sampleFunction
gp_regression_sample(test_x::Union{AbstractArray{Float64, 1}, AbstractArray{Float64, 2}}, n_samples::Int64, gpem::GPModel)

Return n_samples random samples from the Gaussian process posterior, evaluated at test_x. gp_regression.

Arguments

  • test_x: if specified, overrides the test x in gpm. Size $m \times d$.
  • n_samples: integer specifying the number of posterior samples.
  • gpm: the GPModel, that contains the training data (x and y), the kernel, the hyperparameters and the test data for running the regression.
  • full_cov_matrix: whether to use the full covariance matrix or just its diagonal elements (default true).

Return

An array of posterior samples with shape $m \times$ n_samples if n_samples>1 and $m$ otherwise.

source
GpABC.gp_trainMethod
gp_train(gpm::GPModel; <optional keyword arguments>)

Find Maximum Likelihood Estimate of Gaussian Process hyperparameters by maximising gp_loglikelihood, using Optim package. The optimisation target is gp_loglikelihood_log, with gradient computed by gp_loglikelihood_grad. Internally, this function optimises the MLE with respect to logarithms of hyperparameters. This is done for numerical stability. Logarithmisation and exponentiation is performed by this function, i.e. real hyperparameters, not logarithms, are taken in and returned back.

By default, Conjugate Gradient bounded box optimisation is used, as long as the gradient with respect to hyperparameters (covariance_grad) is implemented for the kernel function. If the gradient implementation is not provided, Nelder Mead optimiser is used by default.

Mandatory argument

  • gpm: the GPModel, that contains the training data (x and y), the kernel and the starting hyperparameters that will be used for optimisation.

Optional keyword arguments

  • optimiser::Type{<:Optim.AbstractOptimizer}: the solver to use. If not given, then ConjugateGradient will be used for kernels that have gradient implementation, and NelderMead will be used for those that don't.
  • hp_lower::AbstractArray{Float64, 1}: the lower boundary for box optimisation. Defaults to $e^{-10}$ for all hyperparameters.
  • hp_upper::AbstractArray{Float64, 1}: the upper boundary for box optimisation. Defaults to $e^{10}$ for all hyperparameters.
  • log_level::Int: log level. Default is 0, which is no logging at all. 1 makes gp_train print basic information to standard output. 2 switches Optim logging on, in addition to 1.

Return

The list of all hyperparameters, including the standard deviation of the measurement noise $\sigma_n$. Note that after this function returns, the hyperparameters of gpm will be set to the optimised value, and there is no need to call set_hyperparameters once again.

source