UnivariateApprox

Documentation for UnivariateApprox.

UnivariateApprox.MollifiedBasisType
MollifiedBasis(start, basis, moll)

Take a basis and "mollify" it by a mollifier moll. Only affects basis functions of degree start or higher.

One example is Hermite Functions

Example

julia> basis = ProbabilistHermitePolynomial();

julia> moll = SquaredExponential();

julia> mollified_basis = MollifiedBasis(2, basis, moll); # starts mollifying at quadratic term
source
UnivariateApprox.CreateQuadratureWeightsFunction
CreateQuadratureWeights(pts, functions, true_integrals)

Simple method to create quadrature weights to integrate functions correctly.

Arguments

  • pts::AbstractVector– vector of points
  • functions(N, pt)– A function to evaluate the first N functions at 1d pt
  • true_integrals(N)– A function to return the true integrals of rhe first N functions
  • verbose::Bool– whether to be verbose
source
UnivariateApprox.EvalDiff!Function
EvalDiff!(eval_space::AbstractMatrix, diff_space::AbstractMatrix, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis and its derivative at x and store the results in eval_space and diff_space, respectively.

Example

julia> eval_space = zeros(3,2);

julia> diff_space = zeros(3,2);

julia> EvalDiff!(eval_space, diff_space, LegendrePolynomial(), [0.5, 0.75])

julia> eval_space
3×2 Matrix{Float64}:
  1.0    1.0
  0.5    0.75
 -0.125  0.34375

julia> diff_space
3×2 Matrix{Float64}:
 0.0  0.0
 1.0  1.0
 1.5  2.25

See also: EvalDiff

source
UnivariateApprox.EvalDiffMethod
EvalDiff(max_degree::Int, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis and its derivative at x and return the result.

Example

julia> eval_space, diff_space = EvalDiff(2, LegendrePolynomial(), [0.5, 0.75]);

julia> eval_space
3×2 Matrix{Float64}:
  1.0    1.0
  0.5    0.75
 -0.125  0.34375

julia> diff_space
3×2 Matrix{Float64}:
 0.0  0.0
 1.0  1.0
 1.5  2.25

See also: EvalDiff!

source
UnivariateApprox.EvalDiff2!Function
EvalDiff2!(eval_space::AbstractMatrix, diff_space::AbstractMatrix, diff2_space::AbstractMatrix, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis and its first two derivatives at x and store the results in eval_space, diff_space and diff2_space, respectively.

Example

julia> eval_space = zeros(3,2);

julia> diff_space = zeros(3,2);

julia> diff2_space = zeros(3,2);

julia> EvalDiff2!(eval_space, diff_space, diff2_space, LegendrePolynomial(), [0.5, 0.75])

julia> eval_space
3×2 Matrix{Float64}:
  1.0    1.0
  0.5    0.75
 -0.125  0.34375

julia> diff_space
3×2 Matrix{Float64}:
 0.0  0.0
 1.0  1.0
 1.5  2.25

julia> diff2_space
3×2 Matrix{Float64}:
 0.0  0.0
 0.0  0.0
 3.0  3.0

See also: EvalDiff2

source
UnivariateApprox.EvalDiff2Method
EvalDiff2(max_degree::Int, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis and its first two derivatives at x and return the results.

Example

julia> eval_space, diff_space, diff2_space = EvalDiff2(2, LegendrePolynomial(), [0.5, 0.75]);

julia> eval_space
3×2 Matrix{Float64}:
  1.0    1.0
  0.5    0.75
 -0.125  0.34375

julia> diff_space
3×2 Matrix{Float64}:
 0.0  0.0
 1.0  1.0
 1.5  2.25

julia> diff2_space
3×2 Matrix{Float64}:
 0.0  0.0
 0.0  0.0
 3.0  3.0

See also: EvalDiff2!

source
UnivariateApprox.Evaluate!Function
Evaluate!(space::AbstractMatrix, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis at x and store the result in space.

Example

julia> space = zeros(3,2);

julia> Evaluate!(space, LegendrePolynomial(), [0.5, 0.75])

julia> space
3×2 Matrix{Float64}:
  1.0    1.0
  0.5    0.75
 -0.125  0.34375

See also: Evaluate

source
UnivariateApprox.EvaluateMethod
Evaluate(max_degree::Int, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis at x and return the result.

Example

julia> Evaluate(2, LegendrePolynomial(), [0.5, 0.75])
3×2 Matrix{Float64}:
  1.0    1.0
  0.5    0.75
 -0.125  0.34375

See also: Evaluate!

source
UnivariateApprox.EvaluateDegree!Method
EvaluateDegree!(space::AbstractVector, degree::Int, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis of exactly degree at x and store the result in space.

source
UnivariateApprox.EvaluateDegreeMethod
EvaluateDegree(degree::Int, basis::UnivariateBasis, x::AbstractVector)

Evaluate the univariate basis basis of exactly degree at x and return the result.

source
UnivariateApprox.exactnessFunction
exactness(rule, n)

Return the exactness of the quadrature rule rule for argument n.

Some rules, e.g., Clenshaw-Curtis, have "approximate" exactness beyond the analytical exactness.

source
UnivariateApprox.gausslaguerre_lobattoMethod
gausslaguerre_lobatto(n)

Create an n-point quadrature rule to estimate the integral

\[\int_0^\infty f(t)\exp(-t) dt\]

This method guarantees the first evaluation will be at t_1=0 and is exact on polynomials up to degree 2n-2 (one degree of freedom is lost to fixing t_1).

source
UnivariateApprox.gausspatterson01Method
gausspatterson01(n)

Nested gausspatterson rule to integrate over [0,1], adapted from John Burkardt's implementation.

n must be 1, 3, 7, 15, 31, 63, 127, 255 or 511. Exact for polynomials up to degree 3*2^index - 1 for index > 1.

See here for more information

source