API
Covariance Structures
GaussianRandomFields.CovarianceStructure
— TypeAbstract type CovarianceStructure
Examples
julia> Exponential{Float64} <: CovarianceStructure{Float64}
true
See also: IsotropicCovarianceStructure
, AnisotropicCovarianceStructure
Isotropic Covariance Structures
GaussianRandomFields.IsotropicCovarianceStructure
— TypeAbstract type IsotropicCovarianceStructure <: CovarianceStructure
Examples
julia> Exponential{Float64} <: IsotropicCovarianceStructure{Float64}
true
julia> AnisotropicExponential{Float64} <: IsotropicCovarianceStructure{Float64}
false
See also: Exponential
, Linear
, Spherical
, Whittle
, Gaussian
, SquaredExponential
, Matern
GaussianRandomFields.Linear
— TypeLinear(λ, [σ = 1], [p = 2])
Linear covariance structure with correlation length λ
, (optional) marginal standard deviation σ
and (optional) p
-norm, defined as
$C(x, y) = \begin{cases} σ^2 \left(1 - \displaystyle\frac{ρ}{λ}\right) & \text{if }ρ ≤ λ\\ 0 & \text{if }ρ>λ\end{cases}$
with $ρ = ||x - y||_p$.
Examples
julia> Linear(0.1)
linear (λ=0.1, σ=1.0, p=2.0)
julia> Linear(1.0, σ=2)
linear (λ=1.0, σ=2.0, p=2.0)
See also: Exponential
, Spherical
, Whittle
, Gaussian
, SquaredExponential
, Matern
GaussianRandomFields.Spherical
— TypeSpherical(λ, [σ = 1], [p = 2])
Spherical covariance structure with correlation length λ
, (optional) marginal standard deviation σ
and (optional) p
-norm, defined as
$C(x, y) = \begin{cases} σ^2 \left(1 - \displaystyle\frac{3}{2}\frac{ρ}{λ} + \frac{1}{2}\left(\frac{ρ}{λ}\right)^3\right) & \text{for }ρ≤λ\\0 & \text{for }ρ>λ\end{cases}$
with $ρ = ||x - y||_p$.
Examples
julia> Spherical(0.1)
spherical (λ=0.1, σ=1.0, p=2.0)
julia> Spherical(1.0, σ=2)
spherical (λ=1.0, σ=2.0, p=2.0)
See also: Exponential
, Linear
, Whittle
, Gaussian
, SquaredExponential
, Matern
GaussianRandomFields.Exponential
— TypeExponential(λ, [σ = 1], [p = 2])
Exponential covariance structure with correlation length λ
, (optional) marginal standard deviation σ
and (optional) p
-norm defined as
$C(x, y) = σ^2 \exp\left(-\displaystyle\frac{ρ}{λ}\right)$
with $ρ = ||x - y||_p$.
Examples
julia> Exponential(0.1)
exponential (λ=0.1, σ=1.0, p=2.0)
julia> Exponential(1.0, σ=2)
exponential (λ=1.0, σ=2.0, p=2.0)
See also: Linear
, Spherical
, Whittle
, Gaussian
, SquaredExponential
, Matern
GaussianRandomFields.Whittle
— TypeWhittle(λ, [σ = 1], [p = 2])
Whittle covariance structure with correlation length λ
, (optional) marginal standard deviation σ
and (optional) p
-norm, defined as
$C(x, y) = σ^2 \displaystyle\frac{ρ}{λ} K₁\left(\frac{ρ}{λ}\right)$
with $ρ = ||x-y||_p$.
Examples
julia> Whittle(0.1)
Whittle (λ=0.1, σ=1.0, p=2.0)
julia> Whittle(1.0, σ=2)
Whittle (λ=1.0, σ=2.0, p=2.0)
See also: Exponential
, Linear
, Spherical
, Gaussian
, SquaredExponential
, Matern
GaussianRandomFields.SquaredExponential
— TypeSquaredExponential(λ, [σ = 1], [p = 2])
Squared exponential (Gaussian) covariance structure with correlation length λ
, (optional) marginal standard deviation σ
and (optional) p
-norm, defined as
$C(x, y) = σ^2 \exp\left(-\left(\displaystyle\frac{ρ}{λ}\right)^2\right)$
with $ρ = ||x - y||_p$.
Examples
julia> SquaredExponential(0.1)
Gaussian (λ=0.1, σ=1.0, p=2.0)
julia> SquaredExponential(1, σ=2.)
Gaussian (λ=1.0, σ=2.0, p=2.0)
See also: Exponential
, Linear
, Spherical
, Whittle
, Gaussian
, Matern
GaussianRandomFields.Gaussian
— TypeGaussian(λ, [σ = 1], [p = 2])
Gaussian (squared exponential) covariance structure with correlation length λ
, (optional) marginal standard deviation σ
and (optional) p
-norm, defined as
$C(x, y) = σ \exp\left(-\left(\displaystyle\frac{ρ}{λ}\right)^2\right)$
with $ρ = ||x - y||_p$.
Examples
julia> Gaussian(0.1)
Gaussian (λ=0.1, σ=1.0, p=2.0)
julia> Gaussian(1, σ=2.)
Gaussian (λ=1.0, σ=2.0, p=2.0)
See also: Exponential
, Linear
, Spherical
, Whittle
, SquaredExponential
, Matern
GaussianRandomFields.Matern
— TypeMatern(λ, ν, [σ = 1], [p = 2])
Matérn covariance structure with correlation length λ
, smoothness ν
, (optional) marginal standard deviation σ
and (optional) p
-norm, defined as
$C(x, y) = σ^2 \displaystyle\frac{2^{1 - ν}}{Γ(ν)} \left(\frac{ρ}{λ}\right)^ν K_ν\left(\frac{ρ}{λ}\right)$
with $ρ = ||x - y||_p$.
Examples
julia> Matern(0.1, 1.0)
Matérn (λ=0.1, ν=1.0, σ=1.0, p=2.0)
julia> Matern(1, 1, σ=2.0)
Matérn (λ=1.0, ν=1.0, σ=2.0, p=2.0)
See also: Exponential
, Linear
, Spherical
, Whittle
, Gaussian
, SquaredExponential
Anisotropic Covariance Structures
GaussianRandomFields.AnisotropicCovarianceStructure
— TypeAbstract type AnisotropicCovarianceStructure <: CovarianceStructure
Examples
julia> AnisotropicExponential{Float64} <: AnisotropicCovarianceStructure{Float64}
true
julia> Exponential{Float64} <: AnisotropicCovarianceStructure{Float64}
false
See also: AnisotropicExponential
GaussianRandomFields.AnisotropicExponential
— TypeAnisotropicExponential(A, [σ = 1])
Anisotropic exponential covariance structure with anisotropy matrix A
and (optional) marginal standard deviation σ
, defined as
$C(x, y) = \exp(-ρᵀ A ρ)$
where $ρ = x - y$.
Examples
julia> A = [1 0.5; 0.5 1];
julia> AnisotropicExponential(A)
anisotropic exponential (A=[1.0 0.5; 0.5 1.0], σ=1.0)
Covariance Functions
GaussianRandomFields.AbstractCovarianceFunction
— TypeAbstract type AbstractCovariancFunction
See also: CovarianceFunction
, SeparableCovarianceFunction
GaussianRandomFields.CovarianceFunction
— TypeCovarianceFunction(d, cov)
Covariance function in d
dimensions with covariance structure cov
.
Examples
julia> CovarianceFunction(1, Exponential(0.1))
1d exponential covariance function (λ=0.1, σ=1.0, p=2.0)
julia> CovarianceFunction(2, Matern(0.1, 1.0))
2d Matérn covariance function (λ=0.1, ν=1.0, σ=1.0, p=2.0)
GaussianRandomFields.SeparableCovarianceFunction
— TypeSeparableCovarianceFunction(cov...)
Separable covariance function in length(cov)
dimensions for the covariance structures cov
. Usefull for defining anisotropic covariance functions, or if the non-seperable KarhunenLoeve
expansion is too expensive.
Examples
julia> SeparableCovarianceFunction(Exponential(0.1), Matern(0.01, 1.0))
2d separable covariance function [ exponential (λ=0.1, σ=1.0, p=2.0), Matérn (λ=0.01, ν=1.0, σ=1.0, p=2.0) ]
See also: CovarianceFunction
, KarhunenLoeve
GaussianRandomFields.apply
— Functionapply(cov, pts...)
Returns the covariance matrix, i.e., the covariance function cov
evaluated in the points x
.
Examples
julia> exponential_covariance = CovarianceFunction(1, Exponential(1))
1d exponential covariance function (λ=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=11)
0.0:0.1:1.0
julia> C = apply(exponential_covariance, pts);
julia> heatmap(C);
julia> whittle_covariance = CovarianceFunction(2, Whittle(1))
2d Whittle covariance function (λ=1.0, σ=1.0, p=2.0)
julia> C = apply(whittle_covariance, pts, pts);
julia> heatmap(C);
See also: CovarianceFunction
, Exponential
, Whittle
Gaussian Random Field Generators
GaussianRandomFields.GaussianRandomFieldGenerator
— TypeAbstract type GaussianRandomFieldGenerator
The following Gaussian random field generators are implemented:
Cholesky
: Cholesky factorization of the covariance matrix, exact but expensive for random fields in dimensiond
> 1Spectral
: spectral (eigenvalue) decomposition of the covariance matrix, exact but expensive for random fields in dimensiond
> 1KarhunenLoeve
: Karhunen-Loève expansion, inexact but very efficient for "smooth" random fields when used with a low truncation dimensionCirculantEmbedding
: circulant embedding method, exact and efficient, but can only be used for random fields on structured grids
See also: Cholesky
, Spectral
, KarhunenLoeve
, CirculantEmbedding
Cholesky Factorization
GaussianRandomFields.Cholesky
— TypeCholesky()
Retuns a GaussianRandomFieldGenerator
based on a Cholesky factorization of the covariance matrix.
Examples
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, Cholesky(), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a Cholesky decomposition
julia> heatmap(grf);
The Cholesky factorization requires the covariance matrix to be symmetric and positive definite. If the covariance matrix is not positive definite, an error will be thrown. Try using the Spectral
method in that case.
See also: Spectral
, KarhunenLoeve
, CirculantEmbedding
Spectral Decomposition
GaussianRandomFields.Spectral
— TypeSpectral()
Returns a GaussianRandomFieldGenerator
based on a spectral (eigenvalue) decomposition of the covariance matrix.
Optional Arguments for GaussianRandomField
n::Integer
: the number of eigenvalues to compute. By default, we compute all eigenvalues.eigensolver::EigenSolver
: which method to use for the eigenvalue decomposition (seeAbstractEigenSolver
). The default isEigenSolver()
.
Examples
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, Spectral(), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a spectral decomposition
julia> heatmap(grf);
Try using the Karhunen-Loève expansion if evaluating the covariance matrix is too expensive.
This is also useful when computing Gaussian random fields on a Finite Element mesh using a truncated KL expansion. Here's an example that computes the first 10 eigenfunctions on an L-shaped domain.
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> p, t = Lshape();
julia> grf = GaussianRandomField(cov, Spectral(), p, t, n=10)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a mesh with 998 points and 1861 elements, using a spectral decomposition
See also: Cholesky
, KarhunenLoeve
, CirculantEmbedding
Karhunen Loève Decomposition
GaussianRandomFields.KarhunenLoeve
— TypeKarhunenLoeve(n)
Returns a GaussianRandomFieldGenerator
using a Karhunen-Loève (KL) expansion with n
terms.
Optional Arguments for GaussianRandomField
quad::QuadratureRule
: quadrature rule used for the integral equation (seeQuadratureRule
), default isEOLE()
.nq::Integer
: number of quadrature points in each dimension, where we requirenq^d > n
. Default isnq = ceil(n^(1/d))
.eigensolver::EigenSolver
: which method to use for the eigenvalue decomposition (seeAbstractEigenSolver
). The default isEigenSolver()
.
Examples
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, KarhunenLoeve(300), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 300 terms
julia> plot_eigenvalues(grf); # plot the eigenvalue decay
julia> plot_eigenfunction(grf, 4); # plots the fourth eigenfunction
If more terms n
are used in the expansion, the approximation becomes better.
julia> for n in [1, 2, 5, 10, 20, 50, 100, 200, 500, 1000]
grf = GaussianRandomField(cov, KarhunenLoeve(n), pts, pts)
@printf "%.8f\n" rel_error(grf)
end
0.69838285
0.45494187
0.23231278
0.10079295
0.02647020
0.00978427
0.00356549
0.00107191
0.00019810
0.00004085
Techniqually, the KL expansion is computed using the Nystrom method. For nonstructured grids, we use a bounding box approach. Try using the Spectral
method if this is not what you want.
To avoid an end effect in the eigenvalue decay, choose nq^d ≥ 5n
.
julia> grf = GaussianRandomField(cov, KarhunenLoeve(300), pts, pts, quad=GaussLegendre())
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 300 terms
julia> grf = GaussianRandomField(cov, KarhunenLoeve(300), pts, pts, nq=40)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 300 terms
See also: Cholesky
, Spectral
, CirculantEmbedding
GaussianRandomFields.rel_error
— Functionrel_error(grf)
Returns the relative error in the Karhunen-Loève approximation of the random field, computed as
$1 - \displaystyle\frac{\sum \theta_j^2}{\sigma^2 \int_D \mathrm{d}x}$.
Only useful for fields defined on a rectangular domain.
Examples
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, KarhunenLoeve(300), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 300 terms
julia> @printf "%.8f" rel_error(grf)
0.00046071
Circulant Embedding
GaussianRandomFields.CirculantEmbedding
— TypeCirculantEmbedding()
Returns a GaussianRandomFieldGenerator
that uses FFTs to compute samples of the Gaussian random field.
Circulant embedding can only be applied if the points are specified on a structured grid.
Optional Arguments for GaussianRandomField
minnpadding::Integer
: minimum amount of padding.measure::Bool
: optimize the FFT to increase the efficiency of thesample
method. Default istrue
.primes::Bool
: the size of the minimum circulant embedding of the covariance matrix can be written as a product of small primes (2, 3, 5 and 7). Default isfalse
.
Examples
julia> cov = CovarianceFunction(2, Matern(.1, 1))
2d Matérn covariance function (λ=0.1, ν=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, CirculantEmbedding(), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.1, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a circulant embedding
julia> contourf(grf);
julia> plot_eigenvalues(grf);
With appropriate ordering, the covariance matrix of a Gaussian random field is a (nested block) Toeplitz matrix. This matrix can be embedded into a larger (nested block) circulant matrix, whose eigenvalues can be rapidly computed using FFT. A difficulty here is that although the covariance matrix is positive semi-definite, its circulant extension in general is not. As a remedy, one can add so-called ghost points outside the domain of interest using the optional flag minpadding
.
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, CirculantEmbedding(), pts, pts);
┌ Warning: 318 negative eigenvalues ≥ -0.5828339433508111 detected, Gaussian random field will be
│ approximated (ignoring all negative eigenvalues); increase padding if possible
└ @ GaussianRandomFields ~/.julia/dev/GaussianRandomFields/src/generators/circulant_embedding.jl:94
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a circulant embedding
julia> grf = GaussianRandomField(cov, CirculantEmbedding(), pts, pts, minpadding=79)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a circulant embedding
See also: Cholesky
, Spectral
, KarhunenLoeve
Gaussian Random Fields
GaussianRandomFields.GaussianRandomField
— TypeGaussianRandomField([mean,] cov, generator, pts...)
GaussianRandomField([mean,] cov, generator, nodes, elements)
Compute a Gaussian random field with mean mean
and covariance structure cov
defined in the points pts
, and computed using the Gaussian random field generator generator
.
Examples
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> pts = pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> mean = fill(π, (51, 51));
julia> grf = GaussianRandomField(mean, cov, Cholesky(), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a Cholesky decomposition
If no mean
is specified, a zero-mean Gaussian random field is assumed.
julia> grf = GaussianRandomField(cov, Cholesky(), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a Cholesky decomposition
The Gaussian random field generator generator
can be Cholesky()
, Spectral()
, KarhunenLoeve(n)
(where n
is the number of terms in the expansion), or CirculantEmbedding()
. The points pts
can be specified as arguments of type AbstractVector
, in which case a tensor (Kronecker) product is assumed, or as a Finite Element mesh with node table nodes
and element table elements
.
julia> grf = GaussianRandomField(cov, KarhunenLoeve(500), pts, pts)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 500 terms
julia> exponential_cov = CovarianceFunction(2, Exponential(.1))
2d exponential covariance function (λ=0.1, σ=1.0, p=2.0)
julia> grf = GaussianRandomField(exponential_cov, CirculantEmbedding(), pts, pts)
Gaussian random field with 2d exponential covariance function (λ=0.1, σ=1.0, p=2.0) on a 51×51 structured grid, using a circulant embedding
Separable Gaussian random fields can be defined using SeparableCovarianceFunction
.
julia> separable_cov = SeparableCovarianceFunction(Exponential(.1), Exponential(.1))
2d separable covariance function [ exponential (λ=0.1, σ=1.0, p=2.0), exponential (λ=0.1, σ=1.0, p=2.0) ]
julia> grf = GaussianRandomField(separable_cov, KarhunenLoeve(500), pts, pts)
Gaussian random field with 2d separable covariance function [ exponential (λ=0.1, σ=1.0, p=2.0), exponential (λ=0.1, σ=1.0, p=2.0) ] on a 51×51 structured grid, using a KL expansion with 500 terms
julia> plot_eigenfunction(grf, 3);
We also offer support for anisotropic random fields.
julia> anisotropic_cov = CovarianceFunction(2, AnisotropicExponential([500 400; 400 500]))
2d anisotropic exponential covariance function (A=[500 400; 400 500], σ=1.0)
julia> grf = GaussianRandomField(anisotropic_cov, CirculantEmbedding() , pts, pts)
Gaussian random field with 2d anisotropic exponential covariance function (A=[500 400; 400 500], σ=1.0) on a 51×51 structured grid, using a circulant embedding
julia> heatmap(grf);
For irregular domains, specify the points as matrices containing the nodes
and elements
of a finite element mesh. To compute the value of the random field at the element centers, use the optional keyword mode="center"
.
julia> nodes, elements = Lshape();
julia> grf = GaussianRandomField(cov, Spectral(), nodes, elements)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a mesh with 998 points and 1861 elements, using a spectral decomposition
Alternativly, simply pass a Matrix{T}
of size N
by d
to compute the random field at un unstructured grid defined by the given set of points.
Samples from the random field can be computed using the sample
function.
julia> sample(grf);
See also: Cholesky
, Spectral
, KarhunenLoeve
, CirculantEmbedding
, sample
StatsBase.sample
— Functionsample(grf)
sample(grf[, xi])
Take a sample from the Gaussian random field grf
using the (optional) normally distributed random numbers xi
. The vectorxi
must have appropriate length.
Examples
julia> cov = CovarianceFunction(2, Whittle(.1))
2d Whittle covariance function (λ=0.1, σ=1.0, p=2.0)
julia> pts = pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, CirculantEmbedding(), pts, pts)
Gaussian random field with 2d Whittle covariance function (λ=0.1, σ=1.0, p=2.0) on a 51×51 structured grid, using a circulant embedding
julia> sample(grf);
julia> sample(grf, xi=randn(randdim(grf)));
See also: GaussianRandomField
, Matern
, CovarianceFunction
, [CirculantEmbedding
]
GaussianRandomFields.randdim
— Functionranddim(grf)
Returns the number of random numbers used to sample from the Gaussian random field grf
.
Examples
julia> cov = CovarianceFunction(2, Whittle(.1))
2d Whittle covariance function (λ=0.1, σ=1.0, p=2.0)
julia> pts = pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, KarhunenLoeve(200), pts, pts)
Gaussian random field with 2d Whittle covariance function (λ=0.1, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 200 terms
julia> randdim(grf)
200
Plotting
Standard plotting functions such as plot
and plot!
(for one-dimensional Gaussian random fields), and heatmap
, surface
, contour
and contourf
(for two-dimensional Gaussian random fields) are implemented. There are also some convenience plotting functions defined:
GaussianRandomFields.plot_eigenvalues
— Functionplot_eigenvalues(grf)
Log-log plot of the eigenvalues of the Gaussian random field. Only available for Gaussian random field generators of type Spectral
, KarhunenLoeve
and CirculantEmbedding
.
Examples
julia> cov = CovarianceFunction(2, Gaussian(.1))
2d Gaussian covariance function (λ=0.1, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, KarhunenLoeve(100), pts, pts)
Gaussian random field with 2d Gaussian covariance function (λ=0.1, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 100 terms
julia> plot_eigenvalues(grf);
See also: plot_eigenfunction
GaussianRandomFields.plot_eigenfunction
— Functionplot_eigenfunction(grf, n)
Contour plot of the n
th eigenfunction of the Gaussian random field. Only available for Gaussian random field generators of type Spectral
and KarhunenLoeve
.
Examples
julia> cov = CovarianceFunction(2, Gaussian(.1))
2d Gaussian covariance function (λ=0.1, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=51)
0.0:0.02:1.0
julia> grf = GaussianRandomField(cov, KarhunenLoeve(100), pts, pts)
Gaussian random field with 2d Gaussian covariance function (λ=0.1, σ=1.0, p=2.0) on a 51×51 structured grid, using a KL expansion with 100 terms
julia> plot_eigenfunction(grf, 6); # 6th eigenfunction
See also: plot_eigenfunction
GaussianRandomFields.plot_covariance_matrix
— Functionplot_covariance_matrix(grf[, pts...])
Evaluate the covariance function of the Gaussian random field `grf` in the (optional) points `pts` and plot the result.
Examples
julia> cov = CovarianceFunction(2, Matern(.3, 1))
2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0)
julia> pts = range(0, stop=1, length=11)
0.0:0.1:1.0
julia> grf = GaussianRandomField(cov, CirculantEmbedding(), pts, pts, minpadding=7)
Gaussian random field with 2d Matérn covariance function (λ=0.3, ν=1.0, σ=1.0, p=2.0) on a 11×11 structured grid, using a circulant embedding
julia> plot_covariance_matrix(grf);
julia> pts = range(0, stop=1, length=21)
0.0:0.05:1.0
julia> plot_covariance_matrix(grf, pts, pts);
Unstructured Meshes
GaussianRandomFields.star
— Functionnodes, elements = star()
Returns the 720 nodes
and 1284 elements
of a star.
GaussianRandomFields.Lshape
— Functionnodes, elements = Lshape()
Returns the 998 nodes
and 1861 elements
of an L-shape.
Utilities
GaussianRandomFields.QuadratureRule
— TypeAbstract type QuadratureRule
The following quadrature rules are implemented:
Midpoint::QuadratureRule
: the midpoint ruleTrapezoidal::QuadratureRule
: the trapezoidal ruleSimpson::QuadratureRule
: Simpson's ruleGaussLegendre::QuadratureRule
: Gauss-Legendre quadrature ruleEOLE::QuadratureRule
: expansion-optimal linear estimation
See also: Midpoint
, Trapezoidal
, Simpson
, GaussLegendre
, EOLE
GaussianRandomFields.Midpoint
— TypeGaussianRandomFields.Trapezoidal
— TypeGaussianRandomFields.Simpson
— TypeGaussianRandomFields.GaussLegendre
— TypeGaussianRandomFields.EOLE
— TypeGaussianRandomFields.AbstractEigenSolver
— TypeAbstract type AbstractEigenSolver
The following eigensolvers are implemented: -EigenSolver<:AbstractEigenSolver
: eigenvalue decomposition using eigen -EigsSolver<:AbstractEigenSolver
: eigenvalue decomposition using eigs
See also: EigenSolver
, EigsSolver
GaussianRandomFields.EigenSolver
— TypeGaussianRandomFields.EigsSolver
— Type