GLOM Functionality
GPLinearODEMaker.nlogL_matrix_workspace — TypenlogL_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 onnlogL_hyperparameters∇nlogL_hyperparameters::Vector: Thenloglgradient atnlogL_hyperparametersβs::Vector{Matrix}: List ofinv(Σ) * dΣ_dθwheredΣ_dθis the partial derivative ofΣw.r.t. eachnlogL_hyperparameters
GPLinearODEMaker.GP_posteriors — MethodGP_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
GPLinearODEMaker.GP_posteriors — MethodGP_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
GPLinearODEMaker.GP_posteriors — MethodGP_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
GPLinearODEMaker.GP_posteriors_from_covariances — MethodGP_posteriors_from_covariances(y_obs, Σ_obs, Σ_obs_raw)Calculate the observed posterior mean for the GP used to calculate Σ_obs and Σ_obs_raw
GPLinearODEMaker.GP_posteriors_from_covariances — MethodGP_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
GPLinearODEMaker.GP_posteriors_from_covariances — MethodGP_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
GPLinearODEMaker._covariance! — Methodcovariance!(Σ::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 _.
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 _.
GPLinearODEMaker.calculate_shared_nlogL_matrices! — Methodcalculate_shared_nlogL_matrices!(workspace, glo, non_zero_hyperparameters)Calculates the quantities shared by the nlogL and ∇nlogL calculations and stores them in the existing workspace
GPLinearODEMaker.calculate_shared_nlogL_matrices — Methodcalculate_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
GPLinearODEMaker.calculate_shared_∇nlogL_matrices! — Methodcalculate_shared_∇nlogL_matrices!(workspace, glo, non_zero_hyperparameters)Calculates the quantities shared by the ∇nlogL and ∇∇nlogL calculations and stores them in the existing workspace
GPLinearODEMaker.calculate_shared_∇nlogL_matrices — Methodcalculate_shared_∇nlogL_matrices(glo, non_zero_hyperparameters; kwargs...)Calculates the quantities shared by the ∇nlogL and ∇∇nlogL calculations
GPLinearODEMaker.coefficient_orders — Methodcoefficient_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
GPLinearODEMaker.covariance! — Methodcovariance!(Σ, 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 eachkernel_hyperparameter. (i.e. dorder=[1, 0, 1] would correspond to differenting once w.r.t. the first and thirdkernel_hyperparameter)symmetric::Bool=false: If you know that the resulting covariance matrix should be symmetric, setting this totruecan reduce redundant calculations
GPLinearODEMaker.covariance — Methodcovariance(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
GPLinearODEMaker.covariance — Methodcovariance(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 currentavalues (GLOM coeffients that describe how to combine the differentiated versions of the latent GP) followed by thekernel_hyperparameters(i.e. lengthscales and periods)dΣdθs_total=Int64[]: Which of thetotal_hyperparametersto differentiate the covariance function w.r.t. (i.e. for aglowith 6avalues and 2 kernel hyperparameters, dΣdθs_total=[4, 8] would correspond to differenting once w.r.t. the fourthavalue and second kernel hyperparameter)
See also: covariance!
GPLinearODEMaker.covariance — Methodcovariance(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!
GPLinearODEMaker.covariance_permutations — Methodcovariance_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.
GPLinearODEMaker.covariance_permutations — Methodcovariance_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.
GPLinearODEMaker.d2nlogLdθ — Methodd2nlogLdθ(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 pointy12::Vector: The derivative of observations w.r.t the both parameters at each time pointα::Vector:inv(Σ) * yα1::Vector:inv(Σ) * y1wherey1is the derivative of observations w.r.t the first parameter at each time point
GPLinearODEMaker.d2nlogLdθ — Methodd2nlogLdθ(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 pointy1::Vector: The derivative of observations w.r.t they-affecting parameter at each time pointα::Vector:inv(Σ) * yα1::Vector:inv(Σ) * y1wherey1is the derivative of observations w.r.t the first parameter at each time pointβ2::Matrix:inv(Σ) * dΣ_dθwheredΣ_dθis the partial derivative of theΣw.r.t. theΣ-affecting hyperparameter
GPLinearODEMaker.d2nlogLdθ — Methodd2nlogLdθ(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θwheredΣ_dθis the partial derivative of theΣw.r.t. the first hyperparameterβ2::Matrix: Same asβ1but for the second hyperparameterβ12::Matrix:inv(Σ_obs) * d2Σ_dθ1dθ2whered2Σ_dθ1dθ2is the partial derivative of the covariance matrix Σ_obs w.r.t. both of the hyperparameters being considered
GPLinearODEMaker.dif_coefficients! — Methoddif_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
GPLinearODEMaker.dif_coefficients! — Methoddif_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
GPLinearODEMaker.dnlogLdθ — MethoddnlogLdθ(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θwheredΣ_dθis the partial derivative of theΣw.r.t. a hyperparameter
GPLinearODEMaker.dnlogLdθ — MethoddnlogLdθ(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
GPLinearODEMaker.get_σ — Methodget_σ(glo, x_samp, total_hyperparameters)Calculate the glo GP (using total_hyperparameters) posterior standard deviation at each x_samp point.
GPLinearODEMaker.get_σ — Methodget_σ(L_obs, Σ_obs_samp, diag_Σ_samp)Calculate the GP posterior standard deviation at each sampled point. Algorithm 2.1 from Rasmussen and Williams
GPLinearODEMaker.include_kernel — Methodinclude_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
GPLinearODEMaker.nlogL — MethodnlogL(Σ, 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)
GPLinearODEMaker.nlogL_GLOM! — MethodnlogL_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.
GPLinearODEMaker.nlogL_GLOM — MethodnlogL_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.
GPLinearODEMaker.prep_parallel_covariance — Methodprep_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.
GPLinearODEMaker.prep_parallel_covariance — Methodprep_parallel_covariance(kernel_name; kwargs...)Make it easy to run the covariance calculations on many processors. Makes sure every worker has access to kernel function.
GPLinearODEMaker.reconstruct_total_hyperparameters — Methodreconstruct_total_hyperparameters(glo, hyperparameters)Reinsert the zero coefficients into the non-zero hyperparameter list if needed
GPLinearODEMaker.Σ_observations — MethodΣ_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.
GPLinearODEMaker.Σ_observations — MethodΣ_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
GPLinearODEMaker.Σ_observations — MethodΣ_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.
GPLinearODEMaker.Σ_observations — MethodΣ_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.
GPLinearODEMaker.∇nlogL — Method∇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 ofinv(Σ) * dΣ_dθwheredΣ_dθis the partial derivative ofΣw.r.t. each hyperparameter
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.
GPLinearODEMaker.∇nlogL_GLOM — Method∇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.
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.
GPLinearODEMaker.∇∇nlogL_GLOM — Method∇∇nlogL_GLOM(glo, total_hyperparameters, Σ_obs, y_obs, α, βs)nlogL Hessian for glo using the non-zero total_hyperparameters to construct the covariance matrices.
GPLinearODEMaker.∇∇nlogL_GLOM — Method∇∇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.