Emulators

CalibrateEmulateSample.Emulators.EmulatorType
struct Emulator{FT<:AbstractFloat}

Structure used to represent a general emulator, independently of the algorithm used.

Fields

  • machine_learning_tool::CalibrateEmulateSample.Emulators.MachineLearningTool: Machine learning tool, defined as a struct of type MachineLearningTool.

  • training_pairs::EnsembleKalmanProcesses.DataContainers.PairedDataContainer{FT} where FT<:AbstractFloat: Normalized, standardized, transformed pairs given the Booleans normalize_inputs, standardize_outputs, retained_svd_frac.

  • input_mean::AbstractVector{FT} where FT<:AbstractFloat: Mean of input; length input_dim.

  • normalize_inputs::Bool: If normalizing: whether to fit models on normalized inputs ((inputs - input_mean) * sqrt_inv_input_cov).

  • normalization::Union{Nothing, LinearAlgebra.UniformScaling{FT}, AbstractMatrix{FT}} where FT<:AbstractFloat: (Linear) normalization transformation; size input_dim × input_dim.

  • standardize_outputs::Bool: Whether to fit models on normalized outputs: outputs / standardize_outputs_factor.

  • standardize_outputs_factors::Union{Nothing, AbstractVector{FT}} where FT<:AbstractFloat: If standardizing: Standardization factors (characteristic values of the problem).

  • decomposition::Union{Nothing, LinearAlgebra.SVD}: The singular value decomposition of obs_noise_cov, such that obs_noise_cov = decomposition.U * Diagonal(decomposition.S) * decomposition.Vt. NB: the SVD may be reduced in dimensions.

  • retained_svd_frac::AbstractFloat: Fraction of singular values kept in decomposition. A value of 1 implies full SVD spectrum information.

source
CalibrateEmulateSample.Emulators.EmulatorMethod
Emulator(
    machine_learning_tool::CalibrateEmulateSample.Emulators.MachineLearningTool,
    input_output_pairs::EnsembleKalmanProcesses.DataContainers.PairedDataContainer{FT<:AbstractFloat};
    obs_noise_cov,
    normalize_inputs,
    standardize_outputs,
    standardize_outputs_factors,
    decorrelate,
    retained_svd_frac
) -> CalibrateEmulateSample.Emulators.Emulator{FT} where Float64<:FT<:AbstractFloat

Positional Arguments

  • machine_learning_tool ::MachineLearningTool,
  • input_output_pairs ::PairedDataContainer

Keyword Arguments

  • obs_noise_cov: A matrix/uniform scaling to provide the observational noise covariance of the data - used for data processing (default nothing),
  • normalize_inputs: Normalize the inputs to be unit Gaussian, in the smallest full-rank space of the data (default true),
  • standardize_outputs: Standardize outputs with by dividing by a vector of provided factors (default false),
  • standardize_outputs_factors: If standardizing, the provided dim_output-length vector of factors,
  • decorrelate: Apply (truncated) SVD to the outputs. Predictions are returned in the decorrelated space, (default true)
  • retained_svd_frac: The cumulative sum of singular values retained after output SVD truncation (default 1.0 - no truncation)
source
GaussianProcesses.predictFunction
predict(
    gp::CalibrateEmulateSample.Emulators.GaussianProcess{CalibrateEmulateSample.Emulators.GPJL},
    new_inputs::AbstractArray{FT<:AbstractFloat, 2}
) -> Tuple{Any, Any}

Predict means and covariances in decorrelated output space using Gaussian process models.

source
predict(
    srfi::CalibrateEmulateSample.Emulators.ScalarRandomFeatureInterface,
    new_inputs::AbstractMatrix;
    multithread
) -> Tuple{Any, Any}

Prediction of data observation (not latent function) at new inputs (passed in as columns in a matrix). That is, we add the observational noise into predictions.

source
predict(
    vrfi::CalibrateEmulateSample.Emulators.VectorRandomFeatureInterface,
    new_inputs::AbstractMatrix
) -> Tuple{Any, Any}

Prediction of data observation (not latent function) at new inputs (passed in as columns in a matrix). That is, we add the observational noise into predictions.

source
predict(
    emulator::CalibrateEmulateSample.Emulators.Emulator{FT<:AbstractFloat},
    new_inputs::AbstractArray{FT<:AbstractFloat, 2};
    transform_to_real,
    vector_rf_unstandardize,
    mlt_kwargs...
) -> Tuple{Any, Any}

Makes a prediction using the emulator on new inputs (each new inputs given as data columns). Default is to predict in the decorrelated space.

source
CalibrateEmulateSample.Emulators.standardizeFunction
standardize(
    outputs::AbstractVecOrMat,
    output_covs::Union{AbstractVector{<:AbstractMatrix}, AbstractVector{<:LinearAlgebra.UniformScaling}},
    factors::AbstractVector
) -> Tuple{Any, Union{AbstractVector{<:AbstractMatrix}, AbstractVector{<:LinearAlgebra.UniformScaling}}}

Standardize with a vector of factors (size equal to output dimension).

source
CalibrateEmulateSample.Emulators.reverse_standardizeFunction
reverse_standardize(
    emulator::CalibrateEmulateSample.Emulators.Emulator{FT<:AbstractFloat},
    outputs::AbstractVecOrMat,
    output_covs::AbstractVecOrMat
) -> Tuple{Any, Any}

Reverse a previous standardization with the stored vector of factors (size equal to output dimension). output_cov is a Vector of covariance matrices, such as is returned by svd_reverse_transform_mean_cov.

source
CalibrateEmulateSample.Emulators.svd_transformFunction
svd_transform(
    data::AbstractArray{FT<:AbstractFloat, 2},
    obs_noise_cov::Union{Nothing, AbstractArray{FT<:AbstractFloat, 2}};
    retained_svd_frac
) -> Tuple{Any, Any}

Apply a singular value decomposition (SVD) to the data

  • data - GP training data/targets; size output_dim × N_samples
  • obs_noise_cov - covariance of observational noise
  • truncate_svd - Project onto this fraction of the largest principal components. Defaults to 1.0 (no truncation).

Returns the transformed data and the decomposition, which is a matrix factorization of type LinearAlgebra.SVD.

Note: If F::SVD is the factorization object, U, S, V and Vt can be obtained via F.U, F.S, F.V and F.Vt, such that A = U * Diagonal(S) * Vt. The singular values in S are sorted in descending order.

source
CalibrateEmulateSample.Emulators.svd_reverse_transform_mean_covFunction
svd_reverse_transform_mean_cov(
    μ::AbstractArray{FT<:AbstractFloat, 2},
    σ2::AbstractVector,
    decomposition::LinearAlgebra.SVD
) -> Tuple{Any, Any}

Transform the mean and covariance back to the original (correlated) coordinate system

  • μ - predicted mean; size output_dim × N_predicted_points.
  • σ2 - predicted variance; size output_dim × N_predicted_points. - predicted covariance; size N_predicted_points array of size output_dim × output_dim.
  • decomposition - SVD decomposition of obs_noise_cov.

Returns the transformed mean (size output_dim × N_predicted_points) and variance. Note that transforming the variance back to the original coordinate system results in non-zero off-diagonal elements, so instead of just returning the elements on the main diagonal (i.e., the variances), we return the full covariance at each point, as a vector of length N_predicted_points, where each element is a matrix of size output_dim × output_dim.

source