From 28202a7a2e5740ac4b0cd9d3e92502cef8f91346 Mon Sep 17 00:00:00 2001 From: PieterjanR Date: Tue, 14 Sep 2021 22:09:22 -0700 Subject: [PATCH] Fix for unstructured grids with Cholsky/Spectral --- Project.toml | 2 +- src/covariance_functions/covariance_functions.jl | 3 +++ src/fem.jl | 2 +- test/test_unstructured.jl | 4 ++++ 4 files changed, 9 insertions(+), 2 deletions(-) diff --git a/Project.toml b/Project.toml index 0b6606d..487d951 100644 --- a/Project.toml +++ b/Project.toml @@ -1,6 +1,6 @@ name = "GaussianRandomFields" uuid = "e4b2fa32-6e09-5554-b718-106ed5adafe9" -version = "2.1.3" +version = "2.1.4" [deps] Arpack = "7d9fca2a-8960-54d3-9f78-7d1dccf2cb97" diff --git a/src/covariance_functions/covariance_functions.jl b/src/covariance_functions/covariance_functions.jl index bc76950..4dfa7bf 100644 --- a/src/covariance_functions/covariance_functions.jl +++ b/src/covariance_functions/covariance_functions.jl @@ -185,6 +185,9 @@ function apply(cov::CovarianceStructure, tx::NTuple{2,AbstractMatrix}, Symmetric(C, :U) end +# evaluate for unstructured grids +apply(cov::CovarianceFunction{d}, x::AbstractMatrix{T}, y::AbstractMatrix{N}) where {d, T, N} = apply(cov.cov, (x,x), (x,x)) + # evaluate for KL eigenfunctions function apply(cov::CovarianceStructure, tx::NTuple{2,AbstractMatrix}, y::Tuple) x = first(tx) # select FE nodes diff --git a/src/fem.jl b/src/fem.jl index f44cdb0..88da821 100644 --- a/src/fem.jl +++ b/src/fem.jl @@ -39,7 +39,7 @@ function GaussianRandomField(mean::Vector{<:Real}, cov::CovarianceFunction{d}, m size(pts, 2) == d || throw(DimensionMismatch("second dimension of points must be equal to $(d)")) method isa CirculantEmbedding && throw(ArgumentError("cannot use circulant embedding with an unstructured grid")) - _GaussianRandomField(mean, cov, method, pts', Matrix{Int}(undef,0,0); kwargs...) + _GaussianRandomField(mean, cov, method, copy(pts'), Matrix{Int}(undef,0,0); kwargs...) end GaussianRandomField(cov::CovarianceFunction, method::GaussianRandomFieldGenerator, diff --git a/test/test_unstructured.jl b/test/test_unstructured.jl index 939597a..d7b69ec 100644 --- a/test/test_unstructured.jl +++ b/test/test_unstructured.jl @@ -53,5 +53,9 @@ grf = GaussianRandomField(cov, KarhunenLoeve(128), pts) z = sample(grf) @test isa(grf,GaussianRandomField{KarhunenLoeve{128}}) + grf = GaussianRandomField(cov, Spectral(), pts) + @test isa(grf,GaussianRandomField{Spectral}) + grf = GaussianRandomField(cov, Cholesky(), pts) + @test isa(grf,GaussianRandomField{Cholesky}) end