Gaussian Processes Reference
GpABC
functions for Gaussian Process Regression. See also Gaussian Processes Overview, Gaussian Processes Examples.
Index
GpABC.GPModel
GpABC.GPModel
GpABC.GPModel
GpABC.GPModel
GpABC.gp_loglikelihood
GpABC.gp_loglikelihood
GpABC.gp_loglikelihood_grad
GpABC.gp_loglikelihood_log
GpABC.gp_regression
GpABC.gp_regression
GpABC.gp_regression_sample
GpABC.gp_train
GpABC.set_hyperparameters
Types and Functions
GpABC.GPModel
— TypeGPModel(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.
GpABC.GPModel
— TypeGPModel(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.
GpABC.GPModel
— TypeGPModel
The main type that is used by most functions within the package.
All data matrices are row-major.
Fields
kernel::AbstractGPKernel
: the kernelgp_training_x::AbstractArray{Float64, 2}
: trainingx
. Size: $n \times d$.gp_training_y::AbstractArray{Float64, 2}
: trainingy
. Size: $n \times 1$.gp_test_x::AbstractArray{Float64, 2}
: testx
. 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 togp_loglikelihood
andgp_loglikelihood_grad
GpABC.GPModel
— MethodGPModel(;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.
GpABC.gp_loglikelihood
— Methodgp_loglikelihood(theta::AbstractArray{Float64, 1}, gpm::GPModel)
GpABC.gp_loglikelihood
— Methodgp_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)\]
GpABC.gp_loglikelihood_grad
— Methodgp_loglikelihood_grad(theta::AbstractArray{Float64, 1}, gpem::GPModel)
Gradient of the log likelihood function (gp_loglikelihood_log
) with respect to logged hyperparameters.
GpABC.gp_loglikelihood_log
— Methodgp_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
.
GpABC.gp_regression
— Methodgp_regression(test_x::Union{AbstractArray{Float64, 1}, AbstractArray{Float64, 2}},
gpem::GPModel; <optional keyword arguments>)
GpABC.gp_regression
— Methodgp_regression(gpm::GPModel; <optional keyword arguments>)
Run the Gaussian Process Regression.
Arguments
gpm
: theGPModel
, 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 ingpm
. Size $m \times d$.log_level::Int
(optional): log level. Default is0
, which is no logging at all.1
makesgp_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 tofalse
(i.e. just the variance).batch_size::Int
(optional): Iffull_covariance_matrix
is set tofalse
, 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 totrue
.
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.
GpABC.gp_regression_sample
— Functiongp_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 ingpm
. Size $m \times d$.n_samples
: integer specifying the number of posterior samples.gpm
: theGPModel
, 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 (defaulttrue
).
Return
An array of posterior samples with shape $m \times$ n_samples
if n_samples
>1 and $m$ otherwise.
GpABC.set_hyperparameters
— Methodset_hyperparameters(gpm::GPModel, hypers::AbstractArray{Float64, 1})
Set the hyperparameters of the GPModel
GpABC.gp_train
— Methodgp_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
: theGPModel
, 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, thenConjugateGradient
will be used for kernels that have gradient implementation, andNelderMead
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 is0
, which is no logging at all.1
makesgp_train
print basic information to standard output.2
switchesOptim
logging on, in addition to1
.
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.