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")
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")
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)$")
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")
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)
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")
Function Reference
SparseGridsKit.CCPoints
— TypeCCPoints(domain::Vector{Real})
Generate Clenshaw–Curtis points on interval [domain[1] domain[2]].
Arguments
domain
: Domain
Returns
CCPoints
struct
SparseGridsKit.CCPoints
— Method(p::CCPoints)(n)
Generate Clenshaw-Curtis points on domain.
Arguments
n
: Number of points
SparseGridsKit.CustomLevel
— Method(level::CustomLevel)(l)
Custom rule
Arguments
l
: Level numberfunction
: Function to generate knots fromn
points
Returns
- `
n
: Number of knots
SparseGridsKit.CustomPoints
— TypeCustomPoints(domain::Vector{Real}, knotsfunction::Function)
Generate custom points on interval [domain[1],domain[2]].
Arguments
domain
: Domainknotsfunction
: Function to generate knots fromn
points
Returns
CustomPoints
struct
SparseGridsKit.CustomPoints
— Method(p::CustomPoints)(n)
Generate custom points on domain.
Arguments
n
: Number of points
SparseGridsKit.Doubling
— Method(level::Doubling)(l)
Doubling rule
Arguments
l
: Level number
Returns
- `
n
: Number of knots
SparseGridsKit.GaussHermitePoints
— TypeGaussHermitePoints()
Generate Gauss–Hermite points.
Arguments
None
Returns
GaussHermitePoints
struct
SparseGridsKit.GaussHermitePoints
— Method(p::GaussHermitePoints)(n)
Generate Gauss-Hermite points on interval (-∞,∞).
Arguments
n
: Number of points
SparseGridsKit.GaussLegendrePoints
— TypeGaussLegendrePoints(domain::Vector{Real})
Generate Gauss–Legendre points on interval [domain[1],domain[2]].
Arguments
domain
: Domain
Returns
GaussLegendrePoints
struct
SparseGridsKit.GaussLegendrePoints
— Method(p::GaussLegendrePoints)(n)
Generate Gauss-Legendre points on domain.
Arguments
n
: Number of points
SparseGridsKit.LejaPoints
— TypeLejaPoints([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
: Domainsymmetric
: 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
SparseGridsKit.LejaPoints
— Method(p::LejaPoints)(n)
Generate Leja points on interval [a,b].
Arguments
n
: Number of points
SparseGridsKit.Level
— Method(level::Level)(n::Integer)
Level to knot number mapping
Arguments
l
: Level number
Returns
- `
n
: Number of knots
SparseGridsKit.Points
— MethodPoints(n::Integer)
Generate knots on interval [a,b] using the specified generator function.
Arguments
n
: Number of points
SparseGridsKit.Tripling
— Method(level::Tripling)(l)
Tripling rule
Arguments
l
: Level number
Returns
- `
n
: Number of knots
SparseGridsKit.TwoStep
— Method(level::TwoStep)(l)
TwoStep rule
Arguments
l
: Level number
Returns
- `
n
: Number of knots
SparseGridsKit.UniformPoints
— TypeUniformPoints(domain::Vector{Real})
Generate uniform points on interval [domain[1],domain[2]].
Arguments
domain
: Domain vector [a,b]
Returns
UniformPoints
struct
SparseGridsKit.UniformPoints
— Method(p::UniformPoints)(n)
Generate uniform points on domain.
Arguments
n
: Number of points