API

BayesianSVD.ARCorrelationMethod
ARCorrelation(t; ρ = 1, metric = Euclidean())

Create an AR(1) correlation matrix of type ARCorrelation <: DependentCorrelation <: Correlation.

See also IdentityCorrelation, ExponentialCorrelation, GaussianCorrelation, and MaternCorrelation.

Arguments

  • t: vector of times

Optional Arguments

  • ρ = 1: correlation parameter
  • metric = Euclidean(): metric used for computing the distance between points. All distances in Distances.jl are supported.

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)
X = reduce(hcat,reshape([[x, y] for x = x, y = y], Nx * Ny))'

Ω = SparseCorrelation(x, ρ = 1, metric = Euclidean())
Ω = SparseCorrelation(x, y, ρ = 1, metric = Euclidean())
Ω = SparseCorrelation(X', ρ = 1, metric = Euclidean())
source
BayesianSVD.DataMethod
Data(Y, ulocs, vlocs, k)
Data(Y, X, ulocs, vlocs, k)

Creates the data class.

See also Pars, Posterior, and SampleSVD.

Arguments

  • Y: data of dimension n × m
  • X: covariate matrix of dimension nm × p
  • ulocs: locations for the U basis functions, corresponds to the locations for the rows of Y
  • vlocs: locations for the V basis functions, corresponds to the locations for the columns of Y
  • k: number of basis functions to keep

Examples

k = 5
ΩU = MaternCorrelation(x, ρ = 3, ν = 3.5, metric = Euclidean())
ΩV = MaternCorrelation(t, ρ = 3, ν = 3.5, metric = Euclidean())
data = Data(Z, x, t, k)
pars = Pars(data, ΩU, ΩV)
source
BayesianSVD.ExponentialCorrelationMethod
ExponentialCorrelation(x; ρ = 1, metric = Euclidean())
ExponentialCorrelation(x, y; ρ = 1, metric = Euclidean())
ExponentialCorrelation(X; ρ = 1, metric = Euclidean())

Create an Exponential correlation matrix of type ExponentialCorrelation <: DependentCorrelation <: Correlation.

See also IdentityCorrelation, GaussianCorrelation, MaternCorrelation, and SparseCorrelation.

Arguments

  • x: vector of locations
  • y: vector of locations for a second dimension
  • X: matrix of all pairwise locations

Optional Arguments

  • ρ = 1: length-scale parameter
  • metric = Euclidean(): metric used for computing the distance between points. All distances in Distances.jl are supported.

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)
X = reduce(hcat,reshape([[x, y] for x = x, y = y], Nx * Ny))'

Ω = ExponentialCorrelation(x, ρ = 1, metric = Euclidean())
Ω = ExponentialCorrelation(x, y, ρ = 1, metric = Euclidean())
Ω = ExponentialCorrelation(X', ρ = 1, metric = Euclidean())
source
BayesianSVD.GaussianCorrelationMethod
GaussianCorrelation(x; ρ = 1, metric = Euclidean())
GaussianCorrelation(x, y; ρ = 1, metric = Euclidean())
GaussianCorrelation(X; ρ = 1, metric = Euclidean())

Create an Gaussian correlation matrix of type GaussianCorrelation <: DependentCorrelation <: Correlation.

See also IdentityCorrelation, ExponentialCorrelation, MaternCorrelation, and SparseCorrelation.

Arguments

  • x: vector of locations
  • y: vector of locations for a second dimension
  • X: matrix of all pairwise locations

Optional Arguments

  • ρ = 1: length-scale parameter
  • metric = Euclidean(): metric used for computing the distance between points. All distances in Distances.jl are supported.

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)
X = reduce(hcat,reshape([[x, y] for x = x, y = y], Nx * Ny))'

Ω = GaussianCorrelation(x, ρ = 1, metric = Euclidean())
Ω = GaussianCorrelation(x, y, ρ = 1, metric = Euclidean())
Ω = GaussianCorrelation(X', ρ = 1, metric = Euclidean())
source
BayesianSVD.IdentityCorrelationMethod
IdentityCorrelation(x)
IdentityCorrelation(x, y)
IdentityCorrelation(X)

Create an Identity correlation matrix of type IdentityCorrelation <: IndependentCorrelation <: Correlation.

See also ExponentialCorrelation, GaussianCorrelation, MaternCorrelation, and SparseCorrelation.

Arguments

  • x: vector of locations
  • y: vector of locations for a second dimension
  • X: matrix of all pairwise locations

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)
X = reduce(hcat,reshape([[x, y] for x = x, y = y], Nx * Ny))'

Ω = IdentityCorrelation(x)
Ω = IdentityCorrelation(x, y)
Ω = IdentityCorrelation(X')
source
BayesianSVD.MaternCorrelationMethod
MaternCorrelation(x; ρ = 1, ν = 3.5 metric = Euclidean())
MaternCorrelation(x, y; ρ = 1, ν = 3.5, metric = Euclidean())
MaternCorrelation(X; ρ = 1, ν = 3.5, metric = Euclidean())

Create an Matern correlation matrix of type MaternCorrelation <: DependentCorrelation <: Correlation.

\[K_{\nu, \rho}(s, s') = \frac{2^{1-\nu}}{\Gamma(\nu)}\left(2\nu \frac{||s-s'||}{\rho}\right)^{\nu}J_{\nu}\left(2\nu \frac{||s-s'||}{\rho}\right),\]

for $s,s' \in \mathcal{S}$, where $\Gamma$ is the gamma function, $J_{\nu}$ is the Bessel function of the second kind, and $\{\rho, \nu\}$ are hyperparameters that describe the length-scale and differentiability, respectively.

See also IdentityCorrelation, ExponentialCorrelation, GaussianCorrelation, and SparseCorrelation.

Arguments

  • x: vector of locations
  • y: vector of locations for a second dimension
  • X: matrix of all pairwise locations

Optional Arguments

  • ρ = 1: length-scale parameter
  • ν = 3.5: smoothness parameter
  • metric = Euclidean(): metric used for computing the distance between points. All distances in Distances.jl are supported.

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)
X = reduce(hcat,reshape([[x, y] for x = x, y = y], Nx * Ny))'

Ω = MaternCorrelation(x, ρ = 1, ν = 3.5, metric = Euclidean())
Ω = MaternCorrelation(x, y, ρ = 1, ν = 3.5, metric = Euclidean())
Ω = MaternCorrelation(X', ρ = 1, ν = 3.5, metric = Euclidean())
source
BayesianSVD.ParsMethod
Pars(data::Data, ΩU::Correlation, ΩV::Correlation)

Creates the parameter class of type typeof(data).'

See also Data, Posterior, and SampleSVD.

Arguments

  • data: data structure of type Data
  • ΩU: data structure of type Correlation
  • ΩV: data structure of type Correlation

Examples

k = 5
ΩU = MaternCorrelation(x, ρ = 3, ν = 3.5, metric = Euclidean())
ΩV = MaternCorrelation(t, ρ = 3, ν = 3.5, metric = Euclidean())
data = Data(Z, x, t, k)
pars = Pars(data, ΩU, ΩV)
source
BayesianSVD.PosteriorType
Posterior

Structure of type posterior with subtypes Identity, Exponential, Gaussian, or Matern. Contains the raw posterior samples and some means and 95% quantiles of parameters. Plotting associated with the structure.

See also Pars, Data, and SampleSVD.

Examples

k = 5
ΩU = MaternCorrelation(x, ρ = 3, ν = 3.5, metric = Euclidean())
ΩV = MaternCorrelation(t, ρ = 3, ν = 3.5, metric = Euclidean())
data = Data(Z, x, t, k)
pars = Pars(data, ΩU, ΩV)

posterior, pars = SampleSVD(data, pars; nits = 1000, burnin = 500)

plot(posterior, x, size = (900, 600), basis = 'U', linewidth = 2, c = [:red :green :purple])
plot(posterior, t, size = (900, 500), basis = 'V', linewidth = 2, c = [:red :green :purple])

# for spatial basis functions, provide x and y
plot(posterior, x, y)
source
BayesianSVD.SparseCorrelationMethod
SparseCorrelation(x; ρ = 1, metric = Euclidean())
SparseCorrelation(x, y; ρ = 1, metric = Euclidean())
SparseCorrelation(X; ρ = 1, metric = Euclidean())

Create an Gaussian correlation matrix of type SparseCorrelation <: DependentCorrelation <: Correlation.

See also IdentityCorrelation, ExponentialCorrelation, GaussianCorrelation, and MaternCorrelation.

Arguments

  • x: vector of locations
  • y: vector of locations for a second dimension
  • X: matrix of all pairwise locations

Optional Arguments

  • ρ = 1: length-scale parameter
  • metric = Euclidean(): metric used for computing the distance between points. All distances in Distances.jl are supported.

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)
X = reduce(hcat,reshape([[x, y] for x = x, y = y], Nx * Ny))'

Ω = SparseCorrelation(x, ρ = 1, metric = Euclidean())
Ω = SparseCorrelation(x, y, ρ = 1, metric = Euclidean())
Ω = SparseCorrelation(X', ρ = 1, metric = Euclidean())
source
BayesianSVD.CreateDesignMatrixMethod
CreateDesignMatrix(n::Int, m::Int, Xt, Xs, Xst; intercept = true)

Creates a design matrix.

See also Data, Posterior, and SampleSVD.

Arguments

  • n::Int: Number of spatial locations
  • m::Int: Number of temporal locations
  • Xt: Time only covariates of dimension m x pₘ
  • Xs: Spatial only covariates of dimension n x pₙ
  • Xst: Space-time only covariates of dimension mn x pₘₙ
  • intercept = true: should a global intercept be included

Returns

  • X: Design matrix of dimension mn x (pₘ + pₙ + pₘₙ) or mn x (1 + pₘ + pₙ + pₘₙ)

Examples

m = 75
n = 50

#### Spatial
Xs = hcat(sin.((2*pi .* x) ./ 10), cos.((2*pi .* x) ./ 10))

#### Temporal
Xt = hcat(Vector(t))

#### Spatio-temporal
Xst = rand(Normal(), n*m, 2)

#### Create design matrix
CreateDesignMatrix(n, m, Xt, Xs, Xst, intercept = true)
source
BayesianSVD.GenerateDataMethod

GenerateData(ΣU::Correlation, ΣV::Correlation, D, k, ϵ; SNR = false) GenerateData(ΣU::Vector{T}, ΣV::Vector{T}, D, k, ϵ; SNR = false) where T <: Correlation

Generate random basis functions and data given a correlation structure for $U$ and $V$.

Arguments

  • ΣU::Correlation
  • ΣU::Correlation
  • D: vector of length k
  • k: number of basis functions
  • ϵ: standard devaiation of the noise, if SNR = true then ϵ is the signal to noise ratio

Optional Arguments

  • SNR = false: Boolean for if ϵ is standard deviation (false) or signal to noise ratio value (true)

Returns

  • U: U basis functions
  • V: V basis functions
  • Y: True smooth surface
  • Z: Noisy "observed" surface

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)

ΣU = MaternCorrelation(x, ρ = 3, ν = 3.5, metric = Euclidean())
ΣV = MaternCorrelation(t, ρ = 3, ν = 3.5, metric = Euclidean())


D = [40, 30, 20, 10, 5]
k = 5

Random.seed!(2)
U, V, Y, Z = GenerateData(ΣU, ΣV, D, k, 0.1) # standard deviation
U, V, Y, Z = GenerateData(ΣU, ΣV, D, k, 2, SNR = true) # signal to noise
source
BayesianSVD.PONMethod
PON(n, k, Σ)

Create an n by k orthonormal matrix with covariance matrix Sigma.

Arguments

  • n: Number of locations
  • k: Number of basis functions
  • Σ: covariance matrix

Examples

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)

ΣU = MaternKernel(x, ρ = 3, ν = 3.5, metric = Euclidean())
ΣV = MaternKernel(t, ρ = 3, ν = 3.5, metric = Euclidean())

k = 5
Φ = PON(n, k, ΣU.K)
Ψ = PON(n, k, ΣV.K)
source
BayesianSVD.SampleSVDMethod
SampleSVD(data::Data, pars::Pars; nits = 10000, burnin = 5000, show_progress = true)

Runs the MCMC sampler for the Bayesian SVD model.

See also [`Pars`](@ref), [`Data`](@ref), [`Posterior`](@ref), and [`SampleSVDGrouped`](@ref).

Arguments

  • data: Data structure of type Identity, Exponential, Gaussian, or Matern
  • pars: Parameter structure of type Identity, Exponential, Gaussian, or Matern

Optional Arguments

  • nits = 10000: Total number of posterior samples to compute
  • burnin = 5000: Number of samples discarded as burnin
  • show_progress = true: Indicator on whether to show a progress bar (true) or not (false).

Examples


Random.seed!(2)

m = 100
n = 100
x = range(-5, 5, n)
t = range(0, 10, m)

ΣU = MaternCorrelation(x, ρ = 3, ν = 3.5, metric = Euclidean())
ΣV = MaternCorrelation(t, ρ = 3, ν = 3.5, metric = Euclidean())


D = [40, 30, 20, 10, 5]
k = 5
ϵ = 2

U, V, Y, Z = GenerateData(ΣU, ΣV, D, k, ϵ, SNR = true)


ΩU = MaternCorrelation(x, ρ = 3, ν = 3.5, metric = Euclidean())
ΩV = MaternCorrelation(t, ρ = 3, ν = 3.5, metric = Euclidean())
data = Data(Z, x, t, k)
pars = Pars(data, ΩU, ΩV)

posterior, pars = SampleSVD(data, pars; nits = 1000, burnin = 500)
source
BayesianSVD.SampleSVDGroupedMethod
SampleSVDGrouped(data::Data, pars::Pars; nits = 10000, burnin = 5000, show_progress = true)

Runs the MCMC sampler for the Bayesian SVD model. Exact same documentation as SampleSVD except now the length-scales for U and V
are restricted to be the same.

See also Pars, Data, Posterior, and SampleSVD.

source
BayesianSVD.hpdMethod
hpd(x; p=0.95)

Computes the highest posterior density interval of x.

Arguments

  • x: vector of data.

Optional Arguments

  • p: value between 0 and 1 for the probability level.

Examples

x = [1, 3, 2, 5, .2, 1,9]
hpd(x)
source