Knots

The knot functions take the form

x,w = knots(n)

where n is the number of knots required. The outputs are:

  • x the vector of knots,
  • and w the vector of quadrature weights.

The quadrature weights are generally normalised to sum to one (i.e. representing a probability measure).

These underpin the sparse grid construction.

Knot Functions Basics

Often, one chooses to use Clenshaw–Curtis points which are available using the function CCPoints().

using SparseGridsKit, Plots
p = plot()
for ii in 1:2:7
    plot!(CCPoints(); n=ii)
end
xlabel!("Clenshaw-Curtis Points")
Example block output

Knot functions generally default to the unit interval $[-1,1]$. They can also be mapped to other domains by passing an interval as an arguement. For example, if we wish to use Clenshaw-Curtis points on a parameter domain $\Gamma=[0,100]$ we can do the following.

a = 0
b = 100
p = plot(CCPoints([a,b]))
xlabel!("Clenshaw-Curtis Points")
Example block output

Level Functions

Level functions allow us to define sequences of sets of knots with different growth rates. The level functions are mappings from a level $l$ to a number of knots $n$. For example, we often use the approximate doubling rule

\[n(l) = 2^{l-1} + 1 \quad \forall l > 1, n(1) = 1.\]

using SparseGridsKit
Doubling()(1), Doubling()(2), Doubling()(3)
(1, 3, 5)

Example level functions are plotted below.

using Plots, LaTeXStrings
plot(1:5, Linear().(1:5),label="Linear")
plot(1:5, TwoStep().(1:5),label="TwoStep")
plot!(1:5, Doubling().(1:5),label="Doubling")
plot!(1:5, Tripling().(1:5),label="Tripling")
xlabel!(L"Level $l$")
ylabel!(L"$n(l)$")
Example block output

Later, we will normally use nested sets of points to achieve sparsity in the sparse grid construction. For Clenshaw-Curtis points, nestedness is achieved using the doubling rule.

using SparseGridsKit, Plots
p = plot()
for ii in 1:5
    plot!(CCPoints(); n=Doubling()(ii))
end
xlabel!("Clenshaw-Curtis Points")
ylabel!("Level")
Example block output

More Knot Functions

Functionality is not restricted to knot functions included this package. For example, one could use the quadrature points provided in the FastGaussQuadrature.jl package. These can be wrapped as a CustomPoints to later be used to construct a sparse grid.

using SparseGridsKit, FastGaussQuadrature, Plots

CustomGaussChebyshevPoints = CustomPoints([-1.0,1.0], n->gausschebyshevt(n))
CustomGaussLegendrePoints = CustomPoints([-1.0,1.0], n->gausslegendre(n))

p = plot()
plot!(CCPoints(); n=5)
plot!(UniformPoints(); n=5)
plot!(CustomGaussChebyshevPoints; n=5)
plot!(CustomGaussLegendrePoints; n=5)
Example block output

Leja points are also avaibile. This currently either uses an discrete search based approach to iteratively construct points according to

\[z^{k+1} = \textrm{argmax}_{z\in[-1,1]} v(z) \prod_{i=1}^{k} \textrm{abs}(z-z^i)\]

where $v(z)$ is the weight function. The default weight is $v(z)= \sqrt(\rho(z))$ for $\rho(z)=0.5$ and nonsymmetric points.

using SparseGridsKit, Plots

p = plot(legend=:bottom)
symmetric = false
v(z) = sqrt(0.5)
plot!(LejaPoints(), label="default")
plot!(LejaPoints([-1,1],symmetric,:classic,v), label="nonsymmetric")
symmetric = true
plot!(LejaPoints([-1,1],symmetric,:classic), label="symmetric")
plot!(LejaPoints([-1,1],symmetric,:precomputed), label="precomputed,symmetric")
Example block output

Function Reference

SparseGridsKit.CCPointsType
CCPoints(domain::Vector{Real})

Generate Clenshaw–Curtis points on interval [domain[1] domain[2]].

Arguments

  • domain : Domain

Returns

  • CCPoints struct
source
SparseGridsKit.CustomLevelMethod
(level::CustomLevel)(l)

Custom rule

Arguments

  • l : Level number
  • function : Function to generate knots from n points

Returns

  • `n : Number of knots
source
SparseGridsKit.CustomPointsType
CustomPoints(domain::Vector{Real}, knotsfunction::Function)

Generate custom points on interval [domain[1],domain[2]].

Arguments

  • domain : Domain
  • knotsfunction : Function to generate knots from n points

Returns

  • CustomPoints struct
source
SparseGridsKit.GaussLegendrePointsType
GaussLegendrePoints(domain::Vector{Real})

Generate Gauss–Legendre points on interval [domain[1],domain[2]].

Arguments

  • domain : Domain

Returns

  • GaussLegendrePoints struct
source
SparseGridsKit.LejaPointsType
LejaPoints([a,b],type,v,symmetric)
LejaPoints() = LejaPoints([-1,1], :precomputed, z->sqrt(0.5), false)
LejaPoints(domain) = LejaPoints(domain, :precomputed, z->sqrt(0.5), false)
LejaPoints(domain,type) = LejaPoints(domain,type, z->sqrt(0.5), false)

Generate Leja points on interval [domain[1],domain[2].

Arguments

  • domain : Domain
  • symmetric : Boolean indicating if the points are symmetric (default is false).
  • type : Type of Leja points (:precomputed or :classic), to use optimisation based or discrete search minimisation.
  • v : Function to compute the weight function (default is sqrt(0.5)).

Returns

  • LejaPoints struct
source
SparseGridsKit.LevelMethod
(level::Level)(n::Integer)

Level to knot number mapping

Arguments

  • l : Level number

Returns

  • `n : Number of knots
source
SparseGridsKit.PointsMethod
Points(n::Integer)

Generate knots on interval [a,b] using the specified generator function.

Arguments

  • n : Number of points
source
SparseGridsKit.UniformPointsType
UniformPoints(domain::Vector{Real})

Generate uniform points on interval [domain[1],domain[2]].

Arguments

  • domain : Domain vector [a,b]

Returns

  • UniformPoints struct
source