ApproxFunBase.ArraySpaceType
ArraySpace(s::Space,dims...)

is used to represent array-valued expansions in a space s. The coefficients are of each entry are interlaced.

For example,

f = Fun(x->[exp(x),sin(x)],-1..1)
space(f) == ArraySpace(Chebyshev(),2)
ApproxFunBase.ConversionMethod
Conversion(fromspace::Space, tospace::Space)

Represent a conversion operator between fromspace and tospace, when available.

See also PartialInverseOperator that might be able to represent the inverse, even if this isn't banded.

ApproxFunBase.DefiniteIntegralMethod
DefiniteIntegral([sp::Space])

Return the operator that integrates a Fun over its domain. If sp is unspecified, it is inferred at runtime from the context.

Examples

julia> f = Fun(x -> 3x^2, Chebyshev());

julia> DefiniteIntegral() * f ≈ 2 # integral of 3x^2 over -1..1
true
ApproxFunBase.DerivativeMethod
Derivative(k)

Return the k-th derivative, acting on an unset space. Spaces will be inferred when applying or manipulating the operator. If k is an Int, this returns a derivative in an univariate space. If k is an AbstractVector{Int}, this returns a partial derivative in a multivariate space.

Examples

julia> Derivative(1) * Fun(x->x^2) ≈ Fun(x->2x)
true

julia> Derivative([0,1]) * Fun((x,y)->x^2+y^2) ≈ Fun((x,y)->2y)
true
ApproxFunBase.DerivativeMethod
Derivative(sp::Space, k::AbstractVector{Int})

Return a partial derivative operator on a multivariate space. For example,

Dx = Derivative(Chebyshev()^2,[1,0]) # ∂/∂x
Dy = Derivative(Chebyshev()^2,[0,1]) # ∂/∂y
Tip

Using a static vector as the second argument would help with type-stability.

Examples

julia> ∂y = Derivative(Chebyshev()^2, [0,1]);

julia> ∂y * Fun((x,y)->x^2 + y^2) ≈ Fun((x,y)->2y)
true
ApproxFunBase.DerivativeMethod
Derivative(sp::Space, k::Int)

Return the k-th derivative operator on the space sp.

Examples

julia> Derivative(Chebyshev(), 2) * Fun(x->x^4) ≈ Fun(x->12x^2)
true
ApproxFunBase.DerivativeMethod
Derivative(sp::Space)

Return the first derivative operator, equivalent to Derivative(sp,1).

Examples

julia> Derivative(Chebyshev()) * Fun(x->x^2) ≈ Fun(x->2x)
true
ApproxFunBase.DerivativeMethod
Derivative()

Return the first derivative on an unset space. Spaces will be inferred when applying or manipulating the operator.

Examples

julia> Derivative() * Fun(x->x^2) ≈ Fun(x->2x)
true
ApproxFunBase.DirichletMethod

Dirichlet(sp,k) is the operator associated with restricting the k-th derivative on the boundary for the space sp.

ApproxFunBase.DirichletMethod

Dirichlet(sp) is the operator associated with restricting the the boundary for the space sp.

ApproxFunBase.EvaluationMethod

Evaluation(sp,x,k) is the functional associated with evaluating the k-th derivative at a point x for the space sp.

ApproxFunBase.EvaluationMethod

Evaluation(sp,x) is the functional associated with evaluating at a point x for the space sp.

ApproxFunBase.FunMethod
Fun(f, d::Domain)

Return Fun(f, Space(d)), that is, it uses the default space for the specified domain.

Examples

julia> f = Fun(x->x^2, 0..1);

julia> f(0.1) ≈ (0.1)^2
true
ApproxFunBase.FunMethod
Fun(f, s::Space)

Return a Fun representing the function, number, or vector f in the space s. If f is vector-valued, it Return a vector-valued analogue of s.

Examples

julia> f = Fun(x->x^2, Chebyshev())
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> f(0.1) == (0.1)^2
true
ApproxFunBase.FunMethod
Fun(f)

Return Fun(f, space) by choosing an appropriate space for the function. For univariate functions, space is chosen to be Chebyshev(), whereas for multivariate functions, it is a tensor product of Chebyshev() spaces.

Examples

julia> f = Fun(x -> x^2)
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> f(0.1) == (0.1)^2
true

julia> f = Fun((x,y) -> x + y);

julia> f(0.1, 0.2) ≈ 0.3
true
ApproxFunBase.FunMethod
Fun(s::Space, coefficients::AbstractVector)

Return a Fun with the specified coefficients in the space s

Examples

julia> f = Fun(Fourier(), [1,1]);

julia> f(0.1) == 1 + sin(0.1)
true

julia> f = Fun(Chebyshev(), [1,1]);

julia> f(0.1) == 1 + 0.1
true
ApproxFunBase.FunMethod
Fun(s::Space)

Return Fun(identity, s)

Examples

julia> x = Fun(Chebyshev())
Fun(Chebyshev(), [0.0, 1.0])

julia> x(0.1)
0.1
ApproxFunBase.FunMethod
Fun()

Return Fun(identity, Chebyshev()), which represents the identity function in -1..1.

Examples

julia> f = Fun(Chebyshev())
Fun(Chebyshev(), [0.0, 1.0])

julia> f(0.1)
0.1
ApproxFunBase.IntegralMethod
Integral(k::Int)

Return the k-th integral operator, acting on an unset space. Spaces will be inferred when applying or manipulating the operator.

ApproxFunBase.IntegralMethod
Integral(sp::Space, k::Int)

Return the k-th integral operator on sp. There is no guarantee on normalization.

ApproxFunBase.IntegralMethod
Integral(sp::Space)

Return the first integral operator, equivalent to Integral(sp,1).

ApproxFunBase.IntegralMethod
Intergral()

Return the first integral operator on an unset space. Spaces will be inferred when applying or manipulating the operator.

ApproxFunBase.LaplacianMethod
Laplacian()

Return the laplacian operator on an unset space. Spaces will be inferred when applying or manipulating the operator.

ApproxFunBase.LowRankFunType
LowRankFun(f, space::TensorSpace)

Return an approximation to a bivariate function in a low-rank form

\[f(x,y) = \sum_i \sigma_i \phi_i(x) \psi_i(y)\]

where $\sigma_i$ represent the highest singular values, and $\phi_i(x)$ and $\psi_i(y)$ are orthogonal bases. The summation is truncated after an acceptable tolerance is reached.

Examples

julia> f = (x,y) -> x^2 * y^3;

julia> L = LowRankFun(f, Chebyshev() ⊗ Chebyshev());

julia> L(0.1, 0.2) ≈ f(0.1, 0.2)
true
ApproxFunBase.MultiplicationMethod

Multiplication(f::Fun,sp::Space) is the operator representing multiplication by f on functions in the space sp.

ApproxFunBase.MultiplicationMethod

Multiplication(f::Fun) is the operator representing multiplication by f on an unset space of functions. Spaces will be inferred when applying or manipulating the operator.

ApproxFunBase.PartialInverseOperatorType
PartialInverseOperator(O::Operator, bandwidths = (0, Infinities.ℵ₀))

Return an approximate estimate for inv(O), such that PartialInverseOperator(O) * O is banded, and is approximately I up to a bandwidth that is one less than the sum of the bandwidths of O and PartialInverseOperator(O).

Note

Only upper-triangular operators are supported as of now.

Examples

julia> C = Conversion(Chebyshev(), Ultraspherical(1));

julia> P = PartialInverseOperator(C); # default bandwidth

julia> P * C
TimesOperator : Chebyshev() → Chebyshev()
 1.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  ⋯
  ⋅   1.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  ⋱
  ⋅    ⋅   1.0  0.0  0.0  0.0  0.0  0.0  0.0  0.0  ⋱
  ⋅    ⋅    ⋅   1.0  0.0  0.0  0.0  0.0  0.0  0.0  ⋱
  ⋅    ⋅    ⋅    ⋅   1.0  0.0  0.0  0.0  0.0  0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅   1.0  0.0  0.0  0.0  0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0  0.0  0.0  0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0  0.0  0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0  0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    ⋅   1.0  ⋱
  ⋮    ⋱    ⋱    ⋱    ⋱    ⋱    ⋱    ⋱    ⋱    ⋱   ⋱

julia> P = PartialInverseOperator(C, (0, 4)); # specify an upper bandwidth

julia> P * C
TimesOperator : Chebyshev() → Chebyshev()
 1.0  0.0  0.0  0.0  0.0  0.0  -0.5    ⋅     ⋅     ⋅   ⋅
  ⋅   1.0  0.0  0.0  0.0  0.0   0.0  -1.0    ⋅     ⋅   ⋅
  ⋅    ⋅   1.0  0.0  0.0  0.0   0.0   0.0  -1.0    ⋅   ⋅
  ⋅    ⋅    ⋅   1.0  0.0  0.0   0.0   0.0   0.0  -1.0  ⋅
  ⋅    ⋅    ⋅    ⋅   1.0  0.0   0.0   0.0   0.0   0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅   1.0   0.0   0.0   0.0   0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅    1.0   0.0   0.0   0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅     ⋅    1.0   0.0   0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅     ⋅     ⋅    1.0   0.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅     ⋅     ⋅     ⋅    1.0  ⋱
  ⋅    ⋅    ⋅    ⋅    ⋅    ⋅     ⋅     ⋅     ⋅     ⋅   ⋱
ApproxFunBase.ProductFunMethod
ProductFun(coeffs::AbstractMatrix{T}, sp::AbstractProductSpace; [tol=100eps(T)], [chopping=false]) where {T<:Number}

Represent a bivariate function f in terms of the coefficient matrix coeffs, where the coefficients are obtained using a bivariate transform of the function f in the basis sp.

Examples

julia> P = ProductFun([0 0; 0 1], Chebyshev() ⊗ Chebyshev()) # corresponds to (x,y) -> x*y
ProductFun on Chebyshev() ⊗ Chebyshev()

julia> P(0.1, 0.2) ≈ 0.1 * 0.2
true
ApproxFunBase.ProductFunMethod
ProductFun(M::AbstractVector{<:Fun{<:UnivariateSpace}}, sp::UnivariateSpace)

Represent a bivariate function f(x,y) in terms of the univariate coefficient functions from M. The function f may be reconstructed as

\[f\left(x,y\right)=\sum_{i}M_{i}\left(x\right)b_{i}\left(y\right),\]

where $b_{i}\left(y\right)$ represents the $i$-th basis function for the space sp.

Examples

julia> P = ProductFun([zeros(Chebyshev()), Fun(Chebyshev())], Chebyshev()); # corresponds to (x,y)->x*y

julia> P(0.1, 0.2) ≈ 0.1 * 0.2
true
ApproxFunBase.RowVectorType
RowVector(vector)

A lazy-view wrapper of an AbstractVector, which turns a length-n vector into a 1×n shaped row vector and represents the transpose of a vector (although unlike transpose, the elements are not transposed recursively).

By convention, a vector can be multiplied by a matrix on its left (A * v) whereas a row vector can be multiplied by a matrix on its right (such that RowVector(v) * A = RowVector(transpose(A) * v)). It differs from a 1×n-sized matrix by the facts that its transpose returns a vector and the inner product RowVector(v1) * v2 returns a scalar, but will otherwise behave similarly.

ApproxFunBase.SegmentType
Segment(a,b)

represents a line segment from a to b. In the case where a and b are real and a < b, then this is is equivalent to an Interval(a,b).

ApproxFunBase.SkewSymmetricEigensystemType
SkewSymmetricEigensystem(L::Operator, B::Operator, QuotientSpaceType::Type = QuotientSpace)

Represent the eigensystem L v = λ v subject to B v = 0, where L is skew-symmetric with respect to the standard L2 inner product given the boundary conditions B.

Note

No tests are performed to assert that the operator L is skew-symmetric, and it's the user's responsibility to ensure that the operators are compliant.

The optional argument QuotientSpaceType specifies the type of space to be used to denote the quotient space in the basis recombination process. In most cases, the default choice of QuotientSpace is a good one. In specific instances where B is rank-deficient (e.g. it contains a column of zeros, which typically happens if one of the basis elements already satiafies the boundary conditions), one may need to choose this to be a PathologicalQuotientSpace.

Note

No checks on the rank of B are carried out, and it's up to the user to specify the correct type.

ApproxFunBase.SpaceType
Space{D<:Domain, R}

Abstract supertype of various spaces in which a Fun may be defined, where R represents the type of the basis functions over the domain. Space maps the Domain to the type R.

For example, we have

  • Chebyshev{Interval{Float64}} <: Space{Interval{Float64},Float64}
  • Laurent{PeriodicSegment{Float64}} <: Space{PeriodicSegment{Float64},ComplexF64}
  • Fourier{Circle{ComplexF64}} <: Space{Circle{ComplexF64},Float64}
Note

For now, Space doesn't contain any information about the coefficients

ApproxFunBase.SpaceMethod
(s::Space)(n::Integer, points...)

Evaluate Fun(s, [zeros(n); 1])(points...) efficiently without allocating the vector of coefficients.

Examples

julia> Chebyshev()(1, 0.5)
0.5
ApproxFunBase.SpaceMethod
(s::Space)(n::Integer)

Return a Fun with the coefficients being a sparse representation of [zeros(n); 1]. The result is primarily meant to be evaluated at a specific point.

For orthogonal polynomial spaces, the result will usually represent the n-th basis function.

Examples

julia> Chebyshev()(2) == Fun(Chebyshev(), [0, 0, 1])
true
ApproxFunBase.SpaceOperatorType
SpaceOperator

SpaceOperator is used to wrap other operators and change the domain/range space. This is purely a wrapper. No checks are performed to assess the compatibility of spaces, and no conversions are performed.

Examples

julia> C = Conversion(Chebyshev(), Ultraspherical(0.5));

julia> S = ApproxFunBase.SpaceOperator(C, Chebyshev(), Legendre());

julia> rangespace(S)
Legendre()

julia> domainspace(S)
Chebyshev()
ApproxFunBase.SymmetricEigensystemType
SymmetricEigensystem(L::Operator, B::Operator, QuotientSpaceType::Type = QuotientSpace)

Represent the eigensystem L v = λ v subject to B v = 0, where L is self-adjoint with respect to the standard L2 inner product given the boundary conditions B.

Note

No tests are performed to assert that the operator L is self-adjoint, and it's the user's responsibility to ensure that the operators are compliant.

The optional argument QuotientSpaceType specifies the type of space to be used to denote the quotient space in the basis recombination process. In most cases, the default choice of QuotientSpace is a good one. In specific instances where B is rank-deficient (e.g. it contains a column of zeros, which typically happens if one of the basis elements already satiafies the boundary conditions), one may need to choose this to be a PathologicalQuotientSpace.

Note

No checks on the rank of B are carried out, and it's up to the user to specify the correct type.

ApproxFunBase.TensorSpaceType
TensorSpace(a::Space,b::Space)

represents a tensor product of two 1D spaces a and b. The coefficients are interlaced in lexigraphical order.

For example, consider

Fourier()*Chebyshev()  # returns TensorSpace(Fourier(),Chebyshev())

This represents functions on [-π,π) x [-1,1], using the Fourier basis for the first argument and Chebyshev basis for the second argument, that is, φ_k(x)T_j(y), where

φ_0(x) = 1,
φ_1(x) = sin x,
φ_2(x) = cos x,
φ_3(x) = sin 2x,
φ_4(x) = cos 2x
…

By Choosing (k,j) appropriately, we obtain a single basis:

φ_0(x)T_0(y) (= 1),
φ_0(x)T_1(y) (= y),
φ_1(x)T_0(y) (= sin x),
φ_0(x)T_2(y), …
ApproxFunBase.Conversion_normalizedspaceFunction
Conversion_normalizedspace(S::Space, direction::Val{:forward} = Val(:forward))

Return Conversion(S, normalizedspace(S)). This may be concretely inferred for orthogonal polynomial spaces.

Conversion_normalizedspace(S::Space, ::Val{:backward})

Return Conversion(normalizedspace(S), S). This may be concretely inferred for orthogonal polynomial spaces.

ApproxFunBase.NeumannMethod

Neumann(sp) is the operator associated with restricting the normal derivative on the boundary for the space sp. At the moment it is implemented as Dirichlet(sp,1).

ApproxFunBase.NeumannMethod

Neumann() is the operator associated with restricting the normal derivative on the boundary.

ApproxFunBase.bandmatrices_eigenMethod
bandmatrices_eigen(S::Union{SymmetricEigensystem, SkewSymmetricEigensystem}, n::Integer)

Recast the symmetric/skew-symmetric eigenvalue problem L v = λ v subject to B v = 0 to the generalized eigenvalue problem SA v = λ SB v, where SA and SB are banded operators, and return the n × n matrix representations of SA and SB. If S isa SymmetricEigensystem, the returned matrices will be Symmetric.

Note

No tests are performed to assert that the system is symmetric/skew-symmetric, and it's the user's responsibility to ensure that the operators are compliant.

ApproxFunBase.bvpMethod
bvp(d::Domain, k) = vcat([ldiffbc(d,i) for i=0:div(k,2)-1],
                    [rdiffbc(d,i) for i=0:div(k,2)-1])
bvp(d) = bvp(d,2)

The conditions for the k-th order boundary value problem. See also ldiffbc, rdiffbc, ivp and periodic.

ApproxFunBase.canonicalspaceMethod
canonicalspace(s::Space)

Return a space that is used as a default to implement missing functionality, e.g., evaluation. Implement a Conversion operator or override coefficients to support this.

Examples

julia> ApproxFunBase.canonicalspace(NormalizedChebyshev())
Chebyshev()
ApproxFunBase.choosedomainspaceMethod
choosedomainspace(S::Operator,rangespace::Space)

Return a space ret so that promotedomainspace(S,ret) has the specified range space.

ApproxFunBase.coefficientsMethod
coefficients(cfs::AbstractVector, fromspace::Space, tospace::Space) -> Vector

Convert coefficients in fromspace to coefficients in tospace

Examples

julia> f = Fun(x->(3x^2-1)/2);

julia> coefficients(f, Chebyshev(), Legendre()) ≈ [0,0,1]
true

julia> g = Fun(x->(3x^2-1)/2, Legendre());

julia> coefficients(f, Chebyshev(), Legendre()) ≈ coefficients(g)
true
ApproxFunBase.coefficientsMethod
coefficients(f::Fun, s::Space) -> Vector

Return the coefficients of f in the space s, which may not be the same as space(f).

Examples

julia> f = Fun(x->(3x^2-1)/2);

julia> coefficients(f, Legendre()) ≈ [0, 0, 1]
true
ApproxFunBase.coefficientsMethod
coefficients(f::Fun) -> Vector

Return the coefficients of f, corresponding to the space space(f).

Examples

julia> f = Fun(x->x^2)
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> coefficients(f)
3-element Vector{Float64}:
 0.5
 0.0
 0.5
ApproxFunBase.domainMethod
domain(f::Fun)

Return the domain that f is defined on.

Examples

julia> f = Fun(x->x^2);

julia> domain(f) == ChebyshevInterval()
true

julia> f = Fun(x->x^2, 0..1);

julia> domain(f) == 0..1
true
ApproxFunBase.domainspaceMethod
domainspace(op::Operator)

Return the domain space of op. That is, op*f will first convert f to a Fun in the space domainspace(op) before applying the operator.

ApproxFunBase.eigsMethod
λ, V = eigs(A::Operator, n::Integer; tolerance::Float64=100eps())

Compute the eigenvalues and eigenvectors of the operator A. This is done in the following way:

  • Truncate A into an n×n matrix A₁.
  • Compute eigenvalues and eigenvectors of A₁.
  • Filter for those eigenvectors of A₁ that are approximately eigenvectors

of A as well. The tolerance argument controls, which eigenvectors of the approximation are kept.

ApproxFunBase.eigsMethod
λ, V = eigs(BC::Operator, A::Operator, n::Integer; tolerance::Float64=100eps())

Compute n eigenvalues and eigenvectors of the operator A, subject to the boundary conditions BC.

Examples

julia> #= We compute the spectrum of the second derivative,
          subject to zero boundary conditions.
          We solve this eigenvalue problem in the Chebyshev basis =#

julia> S = Chebyshev();

julia> D = Derivative(S, 2);

julia> BC = Dirichlet(S);

julia> λ, v = ApproxFun.eigs(BC, D, 100);

julia> λ[1:10] ≈ [-(n*pi/2)^2 for n in 1:10] # compare with the analytical result
true
ApproxFunBase.evaluateMethod
evaluate(coefficients::AbstractVector, sp::Space, x)

Evaluate the expansion at a point x that lies in domain(sp). If x is not in the domain, the returned value will depend on the space, and should not be relied upon. See extrapolate to evaluate a function at a value outside the domain.

ApproxFunBase.extrapolateMethod
extrapolate(f::Fun,x)

Return an extrapolation of f from its domain to x.

Examples

julia> f = Fun(x->x^2)
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> extrapolate(f, 2) # 2 lies outside the domain -1..1
4.0
ApproxFunBase.hasconcreteconversion_canonicalMethod
hasconcreteconversion_canonical(sp::Space, ::Val{:forward})

Return Conversion(sp, canonicalspace(sp)) is known statically to be a ConcreteConversion. Assumed to be false by default.

hasconcreteconversion_canonical(sp::Space, ::Val{:backward})

Return Conversion(canonicalspace(sp), sp) is known statically to be a ConcreteConversion. Assumed to be false by default.

ApproxFunBase.interlace!Method

This function implements the algorithm described in:

P. Jain, "A simple in-place algorithm for in-shuffle," arXiv:0805.1598, 2008.
ApproxFunBase.isconvertibleMethod
isconvertible(a::Space, b::Space)

Test whether coefficients may be converted from a to b through a banded Conversion operator.

ApproxFunBase.itransformMethod
itransform(s::Space,coefficients::AbstractVector)

Transform coefficients back to values. Defaults to using canonicalspace as in transform.

Examples

julia> F = Fun(x->x^2, Chebyshev())
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> itransform(Chebyshev(), coefficients(F)) ≈ values(F)
true

julia> itransform(Chebyshev(), [0.5, 0, 0.5])
3-element Vector{Float64}:
 0.75
 0.0
 0.75
ApproxFunBase.ivpMethod
ivp(d::Domain, k) = [ldiffbc(d,i) for i=0:k-1]
ivp(d) = ivp(d,2)

The conditions for the k-th order initial value problem. See also ldiffbc, bvp and periodic.

ApproxFunBase.ncoefficientsMethod
ncoefficients(f::Fun) -> Integer

Return the number of coefficients of a fun

Examples

julia> f = Fun(x->x^2)
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> ncoefficients(f)
3
ApproxFunBase.normalizedspaceFunction
normalizedspace(S::Space)

Return a normalized space (or a direct sum of such spaces) over domain(S). For polynomial spaces, this corresponds to normalized polynomials over the same domain. For a DirectSumSpace, this returns a direct sum over normalized spaces.

ApproxFunBase.pointsMethod
points(f::Fun)

Return a grid of points that f can be transformed into values and back.

Examples

julia> f = Fun(x->x^2);

julia> chebypts(n) = [cos((2i+1)pi/2n) for i in 0:n-1];

julia> points(f) ≈ chebypts(ncoefficients(f))
true
ApproxFunBase.pointsMethod
points(s::Space,n::Integer)

Return a grid of approximately n points, for which a transform exists from values at the grid to coefficients in the space s.

Examples

julia> chebypts(n) = [cos((2i+1)pi/2n) for i in 0:n-1];

julia> points(Chebyshev(), 4) ≈ chebypts(4)
true
ApproxFunBase.promoterangespaceMethod
promoterangespace(S::Operator,sp::Space)

Return the operator S acting on the same space, but now return functions in the specified range space sp

ApproxFunBase.rangespaceMethod
rangespace(op::Operator)

Return the range space of op. That is, op*f will return a Fun in the space rangespace(op), provided f can be converted to a Fun in domainspace(op).

ApproxFunBase.setdomainMethod
setdomain(f::Fun, d::Domain)

Return f projected onto domain.

Note

The new function may differ from the original one, as the coefficients are left unchanged.

Examples

julia> f = Fun(x->x^2);

julia> domain(f) == ChebyshevInterval()
true

julia> g = setdomain(f, 0..1);

julia> domain(g) == 0..1
true

julia> coefficients(f) == coefficients(g)
true
ApproxFunBase.spaceMethod
space(f::Fun)

Return the space of f.

Examples

julia> f = Fun(x->x^2)
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> space(f)
Chebyshev()
ApproxFunBase.supportsinplacetransformMethod
supportsinplacetransform(s::Space)

Trait that states if an inplace transform is possible for the space s. In general, this is possible if transform(s, v) has the same eltype and size as v. By default this is assumed to be false.

New spaces may choose to extend this if the result is known statically.

ApproxFunBase.tocanonicalFunction
tocanonical(d, x)

Map the point x in d to a point in DomainSets.canonicaldomain(d).

Examples

julia> tocanonical(0..0.5, 0)
-1.0

julia> tocanonical(0..0.5, 0.5)
1.0
ApproxFunBase.transformMethod
transform(s::Space, vals)

Transform values on the grid specified by points(s,length(vals)) to coefficients in the space s. Defaults to coefficients(transform(canonicalspace(space),values),canonicalspace(space),space)

Examples

julia> F = Fun(x -> x^2, Chebyshev());

julia> coefficients(F)
3-element Vector{Float64}:
 0.5
 0.0
 0.5

julia> transform(Chebyshev(), values(F)) ≈ coefficients(F)
true

julia> v = map(F, points(Chebyshev(), 4)); # custom grid

julia> transform(Chebyshev(), v)
4-element Vector{Float64}:
 0.5
 0.0
 0.5
 0.0
BandedMatrices.bandwidthsMethod
bandwidths(op::Operator)

Return the bandwidth of op in the form (l,u), where l ≥ 0 represents the number of subdiagonals and u ≥ 0 represents the number of superdiagonals.

Base.:\Method
\(A::Operator,b;tolerance=tol,maxlength=n)

solves a linear equation, usually differential equation, where A is an operator or array of operators and b is a Fun or array of funs. The result u will approximately satisfy A*u = b.

Base.Broadcast.broadcastMethod
f.(x::AbstractVector, y':::AbstractMatrix)

Fast evaluation of a LowRankFun on a cartesian grid x ⨂ y.

Base.chopMethod
chop(f::Fun[, tol = 10eps()]) -> Fun

Reduce the number of coefficients by dropping the tail that is below the specified tolerance.

Examples

julia> f = Fun(Chebyshev(), [1,2,3,0,0,0])
Fun(Chebyshev(), [1, 2, 3, 0, 0, 0])

julia> chop(f)
Fun(Chebyshev(), [1, 2, 3])
Base.getindexMethod
(op::Operator)[k,j]

Return the kth coefficient of op*Fun([zeros(j-1);1],domainspace(op)).

Base.getindexMethod
(op::Operator)[f::Fun]

Construct the operator op * Multiplication(f), that is, it multiplies on the right by f first. Note that op * f is different: it applies op to f.

Examples

julia> x = Fun()
Fun(Chebyshev(), [0.0, 1.0])

julia> D = Derivative()
ConcreteDerivative : ApproxFunBase.UnsetSpace() → ApproxFunBase.UnsetSpace()

julia> Dx = D[x] # construct the operator y -> d/dx * (x * y)
TimesOperator : ApproxFunBase.UnsetSpace() → ApproxFunBase.UnsetSpace()

julia> twox = Dx * x # Evaluate d/dx * (x * x)
Fun(Ultraspherical(1), [0.0, 1.0])

julia> twox(0.1) ≈ 2 * 0.1
true
Base.onesMethod
ones(d::Space)

Return the Fun that represents the function one on the specified space.

Examples

julia> ones(Chebyshev())
Fun(Chebyshev(), [1.0])
Base.strideMethod
stride(f::Fun)

Return the stride of the coefficients, checked numerically

Base.valuesMethod
values(f::Fun)

Return f evaluated at points(f).

Examples

julia> f = Fun(x->x^2)
Fun(Chebyshev(), [0.5, 0.0, 0.5])

julia> values(f)
3-element Vector{Float64}:
 0.75
 0.0
 0.75

julia> map(x->x^2, points(f)) ≈ values(f)
true
Base.zerosMethod
zeros(d::Space)

Return the Fun that represents the function one on the specified space.

Examples

julia> zeros(Chebyshev())
Fun(Chebyshev(), [0.0])
DomainSets.dimensionMethod
dimension(s::Space)

Return the dimension of s, which is the maximum number of coefficients.

LazyArrays.cacheMethod
cache(op::Operator)

Caches the entries of an operator, to speed up multiplying a Fun by the operator.

LinearAlgebra.qrMethod
qr(A::Operator)

returns a cached QR factorization of the Operator A. The result QR enables solving of linear equations: if u=QR, then u approximately satisfies A*u = b.