Multi-fidelity Modelling
The construction offers support for multi-fidelity modelling. A special type of knot function is included which simply maps a level to the same integer, to then be use to specify the model. This has a weight $w=1.0$ and a corresponding polynomial equal to $1$ for all inputs. In practice, for interpolation there is no parameter to be evaluated, but this set up provides the correct multi-fidelity structure using the existing single fidelity construction.
using SparseGridsKit
FidelityPoints()(3), Fidelity()(3)
(([3], [1.0]), 1)
A simple multi-fidelity model considers a constant function $f(x)=y_1^2 + sin(y_2) + y_3$, subject to errors $10^{-\alpha}$. Functions are specified a $f(\vec{\alpha},\vec{y})$ where the vector $\vec{\alpha}$ controls the fidelity and $\vec{y}$ is simply the parameter as usual.
f(y) = y[1]^2 + sin(y[2]) + y[3]
nfid = 1
nparams = 3
f_fidelities(alpha,y) = f(y) + 10^(-alpha[1])
f_fidelities (generic function with 1 method)
The rule that must be used for a fidelity is Fidelity
. This returns $1$ for all input levels as there is only one model per level. The knots are FidelityPoints
, which returns the input value plus a weight zero. The knot function is treated as a special case in the sparse grid construction.
maxfidelity = 5
rule = [fill(Fidelity(),nfid)..., fill(Doubling(),nparams)...]
knots = [fill(FidelityPoints(),nfid)..., fill(CCPoints(),nparams)...]
4-element Vector{Points}:
FidelityPoints(Integer[1, 9223372036854775807])
CCPoints{Float64}([-1.0, 1.0])
CCPoints{Float64}([-1.0, 1.0])
CCPoints{Float64}([-1.0, 1.0])
The domain for FidelityPoints
is by default from 1
to the maximum integer, resulting in the large number above.
We wrap the function $f$ using multifidelityfunctionwrapper
. This allows the user to use the single fidelity sparse grid construction, with function calls separated as $z=[\vec{\alpha},\vec{y}]$. Note that knots
is passed as a parameter: this is used to identify the splitting for the function.
f_wrapped = multifidelityfunctionwrapper(f_fidelities,knots)
MI = create_smolyak_miset(nfid+nparams,2)
sg = create_sparsegrid(MI; knots=knots, rule=rule)
@show get_grid_points(sg)
@show f_eval = f_wrapped.(get_grid_points(sg))
33-element Vector{Float64}:
0.1
1.1
-0.9
0.8071067811865474
-0.6071067811865475
0.9414709848078965
-0.7414709848078965
1.9414709848078966
0.2585290151921035
-0.05852901519210349
⋮
0.5999999999999999
0.010000000000000002
1.01
-0.99
0.8514709848078965
-0.8314709848078965
1.01
1.01
0.001
Finally, interpolation and integration are possible on the sparse grid. Interpolation currently requires dummy values to give a complete parameter vector $[\vec{\alpha}, \vec{y}]$. In the future this may be simplified.
interpolate_on_sparsegrid(sg, f_eval, [[fill(1,nfid)..., fill(0.0,nparams)...]])
pcl = precompute_lagrange_integrals(5, knots, rule)
f_on_grid = f_wrapped.(get_grid_points(sg))
# Test integrate_on_sparsegrid
integral_result = integrate_on_sparsegrid(sg, f_on_grid)
0.33433333333333376
It is also possible to use the adaptive sparse grid construction. A costfunction
can be supplied to give a representative cost for each fidelity. This acts on a vector $\vec{α}$ of the fidelities. We must supply knots
and rules
explictly so that the algorithm can distinguish between model and input dimensions.
knots = [fill(FidelityPoints(),nfid)..., fill(CCPoints(),nparams)...]
rule = [fill(Fidelity(),nfid)..., fill(Doubling(),nparams)...]
(sg, f_on_Z) = adaptive_sparsegrid(f_wrapped, nfid + nparams; maxpts = 100, proftol=1e-4, rule = rule, knots = knots, θ=1e-4, type=:deltaint, costfunction=α->10^prod(α))
(SparseGrid(4, [[1.0, 9.223372036854776e18], [-1.0, 1.0], [-1.0, 1.0], [-1.0, 1.0]], [1 1 … 3 4; 1 2 … 1 1; 1 1 … 1 1; 1 1 … 1 1], [-1, 1, 0, 0, 1], Real[1.0 0.0 0.0 0.0; 1.0 1.0 0.0 0.0; 1.0 -1.0 0.0 0.0; 4.0 0.0 0.0 0.0], [-0.33333333333333326, 0.16666666666666669, 0.16666666666666669, 1.0], Points[FidelityPoints(Integer[1, 9223372036854775807]), CCPoints{Float64}([-1.0, 1.0]), CCPoints{Float64}([-1.0, 1.0]), CCPoints{Float64}([-1.0, 1.0])], Level[Fidelity(), Doubling(), Doubling(), Doubling()], SparseGridsKit.sparse_grid_data{4}(SparseGridsKit.SparseTerm{4}[SparseGridsKit.SparseTerm{4}(([1, 1], [1, 1], [1, 1], [1, 1])), SparseGridsKit.SparseTerm{4}(([1, 1], [2, 1], [1, 1], [1, 1])), SparseGridsKit.SparseTerm{4}(([1, 1], [2, 2], [1, 1], [1, 1])), SparseGridsKit.SparseTerm{4}(([1, 1], [2, 3], [1, 1], [1, 1])), SparseGridsKit.SparseTerm{4}(([4, 1], [1, 1], [1, 1], [1, 1]))], [(1, 1, 1, 1), (1, 2, 1, 1), (1, 3, 1, 1), (4, 1, 1, 1)], [1, 2, 1, 3, 4], [-1, 1, 1, 1, 1], [Dict([2, 1] => 2, [1, 1] => 1, [3, 1] => 3, [4, 1] => 4), Dict([4, 2] => 6, [1, 1] => 1, [3, 5] => 3, [4, 5] => 1, [2, 1] => 2, [3, 1] => 2, [4, 8] => 9, [4, 1] => 2, [4, 6] => 8, [2, 3] => 3…), Dict([4, 2] => 6, [1, 1] => 1, [3, 5] => 3, [4, 5] => 1, [2, 1] => 2, [3, 1] => 2, [4, 8] => 9, [4, 1] => 2, [4, 6] => 8, [2, 3] => 3…), Dict([4, 2] => 6, [1, 1] => 1, [3, 5] => 3, [4, 5] => 1, [2, 1] => 2, [3, 1] => 2, [4, 8] => 9, [4, 1] => 2, [4, 6] => 8, [2, 3] => 3…)], [[[1.0], [2.0], [3.0], [4.0]], [[0.0], [1.0, 0.0, -1.0], [1.0, 0.7071067811865475, 0.0, -0.7071067811865475, -1.0], [1.0, 0.9238795325112867, 0.7071067811865475, 0.38268343236508984, 0.0, -0.38268343236508984, -0.7071067811865475, -0.9238795325112867, -1.0]], [[0.0], [1.0, 0.0, -1.0], [1.0, 0.7071067811865475, 0.0, -0.7071067811865475, -1.0], [1.0, 0.9238795325112867, 0.7071067811865475, 0.38268343236508984, 0.0, -0.38268343236508984, -0.7071067811865475, -0.9238795325112867, -1.0]], [[0.0], [1.0, 0.0, -1.0], [1.0, 0.7071067811865475, 0.0, -0.7071067811865475, -1.0], [1.0, 0.9238795325112867, 0.7071067811865475, 0.38268343236508984, 0.0, -0.38268343236508984, -0.7071067811865475, -0.9238795325112867, -1.0]]], [[1.0, 2.0, 3.0, 4.0], [0.0, 1.0, -1.0, 0.7071067811865475, -0.7071067811865475, 0.9238795325112867, 0.38268343236508984, -0.38268343236508984, -0.9238795325112867], [0.0, 1.0, -1.0, 0.7071067811865475, -0.7071067811865475, 0.9238795325112867, 0.38268343236508984, -0.38268343236508984, -0.9238795325112867], [0.0, 1.0, -1.0, 0.7071067811865475, -0.7071067811865475, 0.9238795325112867, 0.38268343236508984, -0.38268343236508984, -0.9238795325112867]], [[[1.0], [1.0], [1.0], [1.0]], [[1.0], [0.16666666666666669, 0.6666666666666667, 0.16666666666666669], [0.033333333333333354, 0.26666666666666666, 0.4, 0.26666666666666666, 0.033333333333333354], [0.007936507936507936, 0.07310932460800906, 0.13968253968253969, 0.1808589293602449, 0.19682539682539682, 0.1808589293602449, 0.13968253968253969, 0.07310932460800906, 0.007936507936507936]], [[1.0], [0.16666666666666669, 0.6666666666666667, 0.16666666666666669], [0.033333333333333354, 0.26666666666666666, 0.4, 0.26666666666666666, 0.033333333333333354], [0.007936507936507936, 0.07310932460800906, 0.13968253968253969, 0.1808589293602449, 0.19682539682539682, 0.1808589293602449, 0.13968253968253969, 0.07310932460800906, 0.007936507936507936]], [[1.0], [0.16666666666666669, 0.6666666666666667, 0.16666666666666669], [0.033333333333333354, 0.26666666666666666, 0.4, 0.26666666666666666, 0.033333333333333354], [0.007936507936507936, 0.07310932460800906, 0.13968253968253969, 0.1808589293602449, 0.19682539682539682, 0.1808589293602449, 0.13968253968253969, 0.07310932460800906, 0.007936507936507936]]])), [0.1, 1.1, 1.1, 0.0001])
Function Reference
SparseGridsKit.Fidelity
— Method(f:Fidelity)(l)
Fidelity level to knot number mapping
Arguments
l
: Fidelity level
Returns
1
: Fidelity level always returns 1
SparseGridsKit.FidelityPoints
— Method(p::FidelityPoints)(n)
Generate fidelity points
Arguments
n
: Number of points
SparseGridsKit.FidelityPoints
— MethodFidelityPoints()
Generate a struct for fidelity points
Returns
FidelityPoints
struct
SparseGridsKit.fidelitymap
— Methodfidelitymap(y, nfid, nparam)
Splits a grid point y to y[1:nfid],y[nfid+1:nfid+nparam]
Arguments
y
: Grid pointnfid
: Number of fidelity indicesnparam
: Number of parameters
Returns
yfid
: fidelity indicesyparam
: parameters
SparseGridsKit.fidelitypoints
— Methodfidelitypoints(n)
Generate indexes as Natural numbers 1,2,3,... for use as fidelity levels
Arguments
n
: Fidelity index.
Returns
z
: nw
: 1.0
SparseGridsKit.multifidelityfunctionwrapper
— Methodmultifidelityfunctionwrapper(f, nfid, nparam)
Wrap a function for multifidelity evaluation
Arguments
f
: Function to wrapz
: Grid pointnfid
: Number of fidelity indicesnparam
: Number of parameters