GLOM Functionality

GPLinearODEMaker.nlogL_matrix_workspaceType
nlogL_matrix_workspace{T<:Real}

A structure that holds all of the relevant information for calculating nlogL() derivatives. Used to prevent recalculations during optimization.

Arguments

  • nlogL_hyperparameters::Vector: The current hyperparameters
  • Σ_obs::Cholesky: The covariance matrix based on nlogL_hyperparameters
  • ∇nlogL_hyperparameters::Vector: The nlogl gradient at nlogL_hyperparameters
  • βs::Vector{Matrix}: List of inv(Σ) * dΣ_dθ where dΣ_dθ is the partial derivative of Σ w.r.t. each nlogL_hyperparameters
source
GPLinearODEMaker.GP_posteriorsMethod
GP_posteriors(kernel_func, x_obs, y_obs, x_samp, measurement_noise, kernel_hyperparameters; return_mean_obs=false, kwargs...)

Calculate the posterior mean and std at x_samp, (optionally) posterior mean at x_obs and (optionally) the posterior covariance matrix for the GP described by kernel_func and kernel_hyperparameters

source
GPLinearODEMaker.GP_posteriorsMethod
GP_posteriors(glo, x_samp, total_hyperparameters; return_mean_obs=false, y_obs=glo.y_obs, kwargs...)

Calculate the posterior mean and std at x_samp, (optionally) posterior mean at glo.x_obs and (optionally) the posterior covariance matrix for the GP described by glo and total_hyperparameters

source
GPLinearODEMaker.GP_posteriorsMethod
GP_posteriors(glo, total_hyperparameters; y_obs=glo.y_obs, kwargs...)

Calculate the posterior mean at glo.x_obs for the GP described by glo and total_hyperparameters

source
GPLinearODEMaker.GP_posteriors_from_covariancesMethod
GP_posteriors_from_covariances(y_obs, Σ_samp, Σ_obs, Σ_samp_obs, Σ_obs_samp, Σ_obs_raw; return_Σ=true, kwargs...)

Calculate the sampled posterior mean and std, observed posterior mean, and (optionally) the posterior covariance matrix for the GP used to calculate Σ_samp, Σ_obs, Σ_samp_obs, Σ_obs_samp, and Σ_obs_raw

source
GPLinearODEMaker.GP_posteriors_from_covariancesMethod
GP_posteriors_from_covariances(y_obs, Σ_samp, Σ_obs, Σ_samp_obs, Σ_obs_samp; return_Σ=true, kwargs...)

Calculate the sampled posterior mean and std and (optionally) the posterior covariance matrix for the GP used to calculate Σ_samp, Σ_obs, Σ_samp_obs, and Σ_obs_samp

source
GPLinearODEMaker._covariance!Method

covariance!(Σ::Matrix, kernelfunc, x1list, x2list, kernelhyperparameters, dorder, symmetric, samex, equal_spacing)

Fills Σ serially with a covariance matrix by evaluating kernel_func with kernel_hyperparametersfor each pair of x1list and x2list entries. Be careful using a function that starts with a _.

source
GPLinearODEMaker._covariance!Method
_covariance!(Σ::SharedArray, kernel_func, x1list, x2list, kernel_hyperparameters, dorder, symmetric, same_x, equal_spacing)

Fills Σ in parallel with a covariance matrix by evaluating kernel_func with kernel_hyperparametersfor each pair of x1list and x2list entries. Be careful using a function that starts with a _.

source
GPLinearODEMaker.calculate_shared_nlogL_matricesMethod
calculate_shared_nlogL_matrices(glo, non_zero_hyperparameters; Σ_obs=Σ_observations(glo, reconstruct_total_hyperparameters(glo, non_zero_hyperparameters); ignore_asymmetry=true))

Calculates the quantities shared by the nlogL and ∇nlogL calculations

source
GPLinearODEMaker.coefficient_ordersMethod
coefficient_orders(n_out, n_dif; a=ones(n_out, n_dif))

Find the powers that each GLOM coefficient is taken to for each part of the matrix construction before differentiating by any hyperparameters.

Outputs

  • coeff_orders::Array{Int, 6}: filled with integers for what power each coefficient is taken to in the construction of a given block of the total covariance matrix. For example, coeff_orders[1,1,2,3,:,:] would tell you the powers of each coefficient (:,:) that are multiplied by the covariance matrix constructed by evaluating the partial derivative of the kernel (once by t1 and twice by t2) at every pair of time points (2,3) in the construction of the first block of the total covariance matrix (1,1)
  • coeff_coeffs::Array{Int, 4}: Filled with ones anywhere that coeff_orders indicates that coefficients exists to multiply a given covariance matrix for a given block
source
GPLinearODEMaker.covariance!Method
covariance!(Σ, kernel_func, x1list, x2list, kernel_hyperparameters; dorder=zeros(Int64, 2 + length(kernel_hyperparameters)), symmetric=false)

Fills Σ with a covariance matrix by evaluating kernel_func with kernel_hyperparametersfor each pair of x1list and x2list entries.

Keyword Arguments

  • dorder::Vector{T2}=zeros(Int64, 2 + length(kernel_hyperparameters)): How often to differentiate the covariance function w.r.t each kernel_hyperparameter. (i.e. dorder=[1, 0, 1] would correspond to differenting once w.r.t. the first and third kernel_hyperparameter)
  • symmetric::Bool=false: If you know that the resulting covariance matrix should be symmetric, setting this to true can reduce redundant calculations
source
GPLinearODEMaker.covarianceMethod
covariance(glo, total_hyperparameters; kwargs...)

Calculates the total GLOM covariance matrix by combining the latent covariance matrices implied by glo at each glo.x_obs

See also: covariance

source
GPLinearODEMaker.covarianceMethod
covariance(glo, x1list, x2list, total_hyperparameters; dΣdθs_total=Int64[], kwargs...)

Calculates the total GLOM covariance matrix by combining the latent covariance matrices implied by glo for each pair of x1list and x2list entries.

Notable Arguments

  • total_hyperparameters::Vector: The current a values (GLOM coeffients that describe how to combine the differentiated versions of the latent GP) followed by the kernel_hyperparameters (i.e. lengthscales and periods)
  • dΣdθs_total=Int64[]: Which of the total_hyperparameters to differentiate the covariance function w.r.t. (i.e. for a glo with 6 a values and 2 kernel hyperparameters, dΣdθs_total=[4, 8] would correspond to differenting once w.r.t. the fourth a value and second kernel hyperparameter)

See also: covariance!

source
GPLinearODEMaker.covarianceMethod
covariance(kernel_func, x1list, x2list, kernel_hyperparameters; kwargs...)

Calculates a covariance matrix by evaluating kernel_func with kernel_hyperparametersfor each pair of x1list and x2list entries.

See also: covariance!

source
GPLinearODEMaker.covariance_permutationsMethod
covariance_permutations(kernel_func, x_obs, x_samp, measurement_noise, kernel_hyperparameters; return_both=false)

Calculate all of the different versions of the covariance matrices using kernel_func with kernel_hyperparameters between each of the pairs of x_obs and x_samp and themselves.

source
GPLinearODEMaker.covariance_permutationsMethod
covariance_permutations(glo, x_samp, total_hyperparameters; return_both=false)

Calculate all of the different versions of the covariance matrices using the glo GP with total_hyperparameters between each of the pairs of glo.x_obs and x_samp and themselves.

source
GPLinearODEMaker.d2nlogLdθMethod
d2nlogLdθ(y2, y12, α, α1)

Second partial derivative of nlogL(Σ, y) w.r.t. two parameters that affect y

Arguments

  • y2::Vector: The derivative of observations w.r.t the second parameter at each time point
  • y12::Vector: The derivative of observations w.r.t the both parameters at each time point
  • α::Vector: inv(Σ) * y
  • α1::Vector: inv(Σ) * y1 where y1 is the derivative of observations w.r.t the first parameter at each time point
source
GPLinearODEMaker.d2nlogLdθMethod
d2nlogLdθ(y, y1, α, α1, β2)

Second partial derivative of nlogL(Σ, y) w.r.t. a parameter that affects y and a hyperparameter that affects Σ

Arguments

  • y::Vector: The observations at each time point
  • y1::Vector: The derivative of observations w.r.t the y-affecting parameter at each time point
  • α::Vector: inv(Σ) * y
  • α1::Vector: inv(Σ) * y1 where y1 is the derivative of observations w.r.t the first parameter at each time point
  • β2::Matrix: inv(Σ) * dΣ_dθ where dΣ_dθ is the partial derivative of the Σ w.r.t. the Σ-affecting hyperparameter
source
GPLinearODEMaker.d2nlogLdθMethod
d2nlogLdθ(y, α, β1, β2, β12)

Second partial derivative of nlogL(Σ, y) w.r.t. two hyperparameters that affect Σ

Arguments

  • y::Vector: The observations at each time point
  • α::Vector: inv(Σ) * y
  • β1::Matrix: inv(Σ) * dΣ_dθ where dΣ_dθ is the partial derivative of the Σ w.r.t. the first hyperparameter
  • β2::Matrix: Same as β1 but for the second hyperparameter
  • β12::Matrix: inv(Σ_obs) * d2Σ_dθ1dθ2 where d2Σ_dθ1dθ2 is the partial derivative of the covariance matrix Σ_obs w.r.t. both of the hyperparameters being considered
source
GPLinearODEMaker.dif_coefficients!Method
dif_coefficients!(n_out, n_dif, dΣdθs_total::Vector, coeff_orders, coeff_coeffs)

Modify coeff_orders and coeff_coeffs with the coefficients for constructing differentiated version of the kernel (for the differentiations implied by dΣdθ_totals) using the powers that each coefficient is taken to for each part of the matrix construction

source
GPLinearODEMaker.dif_coefficients!Method
dif_coefficients!(n_out, n_dif, dΣdθ_total::Int, coeff_orders, coeff_coeffs)

Modify coeff_orders and coeff_coeffs with the coefficients for constructing differentiated version of the kernel (for the differentiation implied by dΣdθ_total) using the powers that each coefficient is taken to for each part of the matrix construction

source
GPLinearODEMaker.dnlogLdθMethod
dnlogLdθ(y, α, β)

First partial derivative of nlogL(Σ, y) w.r.t. hyperparameters that affect Σ

Arguments

  • y::Vector: The observations at each time point
  • α::Vector: inv(Σ) * y
  • β::Matrix: inv(Σ) * dΣ_dθ where dΣ_dθ is the partial derivative of the Σ w.r.t. a hyperparameter
source
GPLinearODEMaker.dnlogLdθMethod
dnlogLdθ(y1, α)

First partial derivative of nlogL(Σ, y) w.r.t. hyperparameters that affect y

Arguments

  • y1::Vector: The derivative of observations at each time point
  • α::Vector: inv(Σ) * y
source
GPLinearODEMaker.get_σMethod
get_σ(glo, x_samp, total_hyperparameters)

Calculate the glo GP (using total_hyperparameters) posterior standard deviation at each x_samp point.

source
GPLinearODEMaker.get_σMethod
get_σ(L_obs, Σ_obs_samp, diag_Σ_samp)

Calculate the GP posterior standard deviation at each sampled point. Algorithm 2.1 from Rasmussen and Williams

source
GPLinearODEMaker.include_kernelMethod
include_kernel(kernel_name)

Tries to include the specified kernel, assuming it was included with GLOM. Returns the kernel function and number of hyperparameters it uses

source
GPLinearODEMaker.nlogLMethod
nlogL(Σ, y; α= Σ \ y, nlogL_normalization=logdet(Σ)+length(y)*log(2*π))

Negative log likelihood for data y assuming it was drawn from a multivariate normal distribution with 0 mean and covariance Σ (usually Σ + noise)

source
GPLinearODEMaker.nlogL_GLOM!Method
nlogL_GLOM!(workspace, glo, total_hyperparameters; y_obs=copy(glo.y_obs))

Negative log likelihood for glo using total_hyperparameters to construct the covariance matrix and storing the intermediate results in workspace.

source
GPLinearODEMaker.nlogL_GLOMMethod
nlogL_GLOM(glo, total_hyperparameters; y_obs=copy(glo.y_obs), kwargs...)

Negative log likelihood for glo using the non-zero total_hyperparameters to construct the covariance matrix.

source
GPLinearODEMaker.prep_parallel_covarianceMethod
prep_parallel_covariance(kernel_name, kernel_path; kwargs...)

Make it easy to run the covariance calculations on many processors. Makes sure every worker has access to kernel function, importing it from the kernel_path.

source
GPLinearODEMaker.Σ_observationsMethod
Σ_observations(kernel_func, x_obs, measurement_noise, kernel_hyperparameters; ignore_asymmetry=false, return_both)

Calculates the covariance matrix of kernel_func at x_obs and adds measurement_noise^2 to the diagonal and performs a Cholesky factorization. Optionally returns the Σ back as well.

source
GPLinearODEMaker.Σ_observationsMethod
Σ_observations(glo, total_hyperparameters; ignore_asymmetry=false, return_both=false)

Calculates a Cholesky decomposition of the GLOM covariance matrix implied by glo and total_hyperparameters including glo.noise or glo.covariance on the (block) diagonal

source
GPLinearODEMaker.Σ_observationsMethod
Σ_observations(Σ, measurement_noise::Vector; return_both=false, kwargs...)

Add measurement_noise^2 to the diagonal of Σ and performs a Cholesky factorization. Optionally returns the Σ back as well.

source
GPLinearODEMaker.Σ_observationsMethod
Σ_observations(Σ, measurement_covariance::Array{T, 3}; return_both=false, kwargs...)

Add measurement_covariance to the block diagonal of Σ and performs a Cholesky factorization. Optionally returns the Σ back as well.

source
GPLinearODEMaker.∇nlogLMethod
∇nlogL(y, α, βs)

Gradient of nlogL(Σ, y) w.r.t. hyperparameters that affect Σ

Arguments

  • y::Vector: The observations at each time point
  • α::Vector: inv(Σ) * y
  • βs::Vector{Matrix}: List of inv(Σ) * dΣ_dθ where dΣ_dθ is the partial derivative of Σ w.r.t. each hyperparameter
source
GPLinearODEMaker.∇nlogL_GLOM!Method
∇nlogL_GLOM!(workspace, glo, total_hyperparameters; y_obs=copy(glo.y_obs))

nlogL gradient for glo using the non-zero total_hyperparameters to construct the covariance matrices and storing the intermediate results in workspace.

source
GPLinearODEMaker.∇nlogL_GLOMMethod
∇nlogL_GLOM(glo, total_hyperparameters; y_obs=copy(glo.y_obs), kwargs...)

nlogL gradient for glo using the non-zero total_hyperparameters to construct the covariance matrices.

source
GPLinearODEMaker.∇∇nlogL_GLOM!Method
∇∇nlogL_GLOM!(workspace, glo, total_hyperparameters; y_obs=copy(glo.y_obs))

nlogL Hessian for glo using the non-zero total_hyperparameters to construct the covariance matrices and storing the intermediate results in workspace.

source
GPLinearODEMaker.∇∇nlogL_GLOMMethod
∇∇nlogL_GLOM(glo, total_hyperparameters, Σ_obs, y_obs, α, βs)

nlogL Hessian for glo using the non-zero total_hyperparameters to construct the covariance matrices.

source
GPLinearODEMaker.∇∇nlogL_GLOMMethod
∇∇nlogL_GLOM(glo, total_hyperparameters; y_obs=copy(glo.y_obs))

nlogL Hessian for glo using the non-zero total_hyperparameters to construct the covariance matrices.

source