Special operators
In this section, we present several operations common to dynamical systems.
Projection
When working with SequenceSpace
, one frequently needs to adjust the order of truncation of the chosen basis. This operation is implemented as the project
and project!
functions. In fact, these functions provide a general mechanism to retrieve a finite part of the infinite dimensional operators introduced later in this section.
Each project
or project!
call verifies a compatibility criterion between spaces. For Sequence
and LinearOperator
, two VectorSpace
are compatible if:
- all comprised
SequenceSpace
only differ from their order. For instance,Taylor(n)
andTaylor(m)
are compatible for any positiven::Int
andm::Int
. However,Taylor(n)
andTensorSpace(Taylor(m), Fourier(k, 1.0))
are not compatible for any positiven::Int
,m::Int
andk::Int
. - all comprised
CartesianSpace
have the same number of cartesian products. For instance,CartesianPower(a, 2)
andCartesianProduct(a, a)
are compatible for anya::VectorSpace
. However,CartesianProduct(a, b)
andCartesianProduct(CartesianPower(a, 1), b)
are not compatible for anya::VectorSpace
andb::VectorSpace
.
julia> A = LinearOperator(Taylor(1) ⊗ Chebyshev(1), Taylor(1) ⊗ Chebyshev(1), [1 0 0 0 ; 0 1 0 0 ; 0 0 1 0 ; 0 0 0 1]) # project(I, Taylor(1) ⊗ Chebyshev(1), Taylor(1) ⊗ Chebyshev(1))
LinearOperator : Taylor(1) ⊗ Chebyshev(1) → Taylor(1) ⊗ Chebyshev(1) with coefficients Matrix{Int64}: 1 0 0 0 0 1 0 0 0 0 1 0 0 0 0 1
julia> project(A, Taylor(1) ⊗ Chebyshev(2), Taylor(2) ⊗ Chebyshev(1))
LinearOperator : Taylor(1) ⊗ Chebyshev(2) → Taylor(2) ⊗ Chebyshev(1) with coefficients Matrix{Int64}: 1 0 0 0 0 0 0 1 0 0 0 0 0 0 0 0 0 0 0 0 1 0 0 0 0 0 0 1 0 0 0 0 0 0 0 0
Moreover, the following identifications are permitted:
julia> a = Sequence(Taylor(1), [1, 1]) # 1 + x
Sequence in Taylor(1) with coefficients Vector{Int64}: 1 1
julia> A = project(a, ParameterSpace(), Taylor(2))
LinearOperator : 𝕂 → Taylor(2) with coefficients Matrix{Int64}: 1 1 0
julia> project(A, space(a))
Sequence in Taylor(1) with coefficients Vector{Int64}: 1 1
Multiplication
Let $V$ be a SequenceSpace
with discrete convolution $*$ and $a \in V$. The multiplication operator Multiplication
represents the mapping $\mathcal{M}_a : V \to V$ defined by
\[\mathcal{M}_a (b) \bydef a * b, \qquad \text{for all } b \in V.\]
The action of Multiplication
is performed by the right product *
of a Multiplication
with a Sequence{<:SequenceSpace}
; alternatively, Multiplication
defines a method on a Sequence{<:SequenceSpace}
representing *
.
julia> a = Sequence(Taylor(1), [1, 1]); # 1 + x
julia> b = Sequence(Taylor(2), [0, 0, 1]); # x^2
julia> a * b
Sequence in Taylor(3) with coefficients Vector{Int64}: 0 0 1 1
julia> ℳ = Multiplication(a)
Multiplication{Sequence{Taylor, Vector{Int64}}}(Sequence(Taylor(1), [1, 1]))
julia> ℳ * b # ℳ(b)
Sequence in Taylor(3) with coefficients Vector{Int64}: 0 0 1 1
A finite dimensional truncation of Multiplication
may be obtained via project
or project!
.
julia> project(ℳ, Taylor(2), image(*, Taylor(1), Taylor(2)))
LinearOperator : Taylor(2) → Taylor(3) with coefficients Matrix{Int64}: 1 0 0 1 1 0 0 1 1 0 0 1
Derivation and integration
Both Derivative
and Integral
have a field order::Union{Int,Tuple{Vararg{Int}}}
to specify how many times the operator is composed with itself. No derivation or integration is performed whenever a value of 0
is given.
julia> a = Sequence(Taylor(2), [1, 1, 1]); # 1 + x + x^2
julia> differentiate(a)
Sequence in Taylor(1) with coefficients Vector{Int64}: 1 2
julia> 𝒟 = Derivative(1)
Derivative{Int64}(1)
julia> 𝒟 * a # 𝒟(a)
Sequence in Taylor(1) with coefficients Vector{Int64}: 1 2
A finite dimensional truncation of Derivative
and Integral
may be obtained via project
or project!
:
julia> project(Derivative(1), Taylor(2), image(Derivative(1), Taylor(2)), Float64)
LinearOperator : Taylor(2) → Taylor(1) with coefficients SparseArrays.SparseMatrixCSC{Float64, Int64}: ⋅ 1.0 ⋅ ⋅ ⋅ 2.0
julia> project(Integral(1), Taylor(2), image(Integral(1), Taylor(2)), Float64)
LinearOperator : Taylor(2) → Taylor(3) with coefficients SparseArrays.SparseMatrixCSC{Float64, Int64}: ⋅ ⋅ ⋅ 1.0 ⋅ ⋅ ⋅ 0.5 ⋅ ⋅ ⋅ 0.3333333333333333
Evaluation
The evaluation operator Evaluation
has a field value::Union{Number,Nothing,Tuple{Vararg{Union{Number,Nothing}}}}
representing the evaluation point. No scaling is performed whenever a value of nothing
is given.
julia> a = Sequence(Taylor(2), [1, 1, 1]); # 1 + x + x^2
julia> evaluate(a, 0.1)
1.11
julia> ℰ = Evaluation(0.1)
Evaluation{Float64}(0.1)
julia> ℰ * a # ℰ(a)
1.11
julia> b = Sequence(Taylor(1) ⊗ Fourier(1, 1.0), [0.5, 0.5, 0.0, 0.0, 0.5, 0.5]); # (1 + x) cos(y)
julia> evaluate(b, (0.1, nothing)) # Evaluation(0.1, nothing) * b
Sequence in Taylor(0) ⊗ Fourier(1, 1.0) with coefficients Vector{Float64}: 0.55 0.0 0.55
Moreover, Evaluation
is defined on CartesianSpace
by acting component-wise.
julia> c = Sequence(Taylor(1)^2, [1, 1, 2, 2]); # 1 + x, 2 + 2x
julia> evaluate(c, 0.1) # Evaluation(0.1) * c
2-element Vector{Float64}: 1.1 2.2
A finite dimensional truncation of Evaluation
may be obtained via project
or project!
:
julia> project(Evaluation(0.1), Taylor(2), image(Evaluation(0.1), Taylor(2)), Float64)
LinearOperator : Taylor(2) → Taylor(0) with coefficients Matrix{Float64}: 1.0 0.1 0.010000000000000002
Furthermore, in the context of Evaluation
, the concept of compatibility between two VectorSpace
is more permissive to allow manipulating Evaluation
more like a functional:
julia> project(Evaluation(0.1), Taylor(2), ParameterSpace(), Float64)
LinearOperator : Taylor(2) → 𝕂 with coefficients Matrix{Float64}: 1.0 0.1 0.010000000000000002
Scale
The scale operator Scale
has a field value::Union{Number,Tuple{Vararg{Number}}}
representing the scaling factor. No scaling is performed whenever a value of 1
is given.
julia> a = Sequence(Taylor(2), [1, 1, 1]) # 1 + x + x^2
Sequence in Taylor(2) with coefficients Vector{Int64}: 1 1 1
julia> scale(a, 2)
Sequence in Taylor(2) with coefficients Vector{Int64}: 1 2 4
julia> 𝒮 = Scale(2)
Scale{Int64}(2)
julia> 𝒮 * a # 𝒮(a)
Sequence in Taylor(2) with coefficients Vector{Int64}: 1 2 4
A finite dimensional truncation of Scale
may be obtained via project
or project!
:
julia> project(Scale(2), Taylor(2), image(Scale(2), Taylor(2)), Float64)
LinearOperator : Taylor(2) → Taylor(2) with coefficients SparseArrays.SparseMatrixCSC{Float64, Int64}: 1.0 ⋅ ⋅ ⋅ 2.0 ⋅ ⋅ ⋅ 4.0
Shift
The shift operator Shift
has a field value::Union{Number,Tuple{Vararg{Number}}}
representing the shift. No shift is performed whenever a value of 0
is given.
Currently, only Fourier
space allows values different than 0
.
julia> a = Sequence(Fourier(1, 1.0), [0.5, 0.0, 0.5]) # cos(x)
Sequence in Fourier(1, 1.0) with coefficients Vector{Float64}: 0.5 0.0 0.5
julia> shift(a, π)
Sequence in Fourier(1, 1.0) with coefficients Vector{ComplexF64}: -0.5 - 6.123233995736766e-17im 0.0 + 0.0im -0.5 + 6.123233995736766e-17im
julia> 𝒮 = Shift(π)
Shift{Irrational{:π}}(π)
julia> 𝒮 * a # 𝒮(a)
Sequence in Fourier(1, 1.0) with coefficients Vector{ComplexF64}: -0.5 - 6.123233995736766e-17im 0.0 + 0.0im -0.5 + 6.123233995736766e-17im
A finite dimensional truncation of Shift
may be obtained via project
or project!
:
julia> project(Shift(π), Fourier(1, 1.0), image(Shift(π), Fourier(1, 1.0)), Complex{Float64})
LinearOperator : Fourier(1, 1.0) → Fourier(1, 1.0) with coefficients SparseArrays.SparseMatrixCSC{ComplexF64, Int64}: -1.0 - 1.2246467991473532e-16im … ⋅ ⋅ ⋅ ⋅ -1.0 + 1.2246467991473532e-16im
API
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, A::LinearOperator)
Represent A
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, J::UniformScaling)
Represent J
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator{ParameterSpace,<:VectorSpace}, a::Sequence)
Represent a
as a LinearOperator
from ParameterSpace
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project
.
RadiiPolynomial.project!
— Methodproject!(c::Sequence, A::LinearOperator{ParameterSpace,<:VectorSpace})
Represent A
as a Sequence
in space(c)
. The result is stored in c
by overwriting it.
See also: project
.
RadiiPolynomial.project!
— Methodproject!(c::Sequence, a::Sequence)
Represent a
as a Sequence
in space(c)
. The result is stored in c
by overwriting it.
See also: project
.
RadiiPolynomial.project
— Methodproject(A::LinearOperator, domain_dest::VectorSpace, codomain_dest::VectorSpace, ::Type{T}=eltype(A))
Represent A
as a LinearOperator
from domain_dest
to codomain_dest
.
See also: project!
.
RadiiPolynomial.project
— Methodproject(A::LinearOperator{ParameterSpace,<:VectorSpace}, space_dest::VectorSpace, ::Type{T}=eltype(A))
Represent A
as a Sequence
in space_dest
.
See also: project!
.
RadiiPolynomial.project
— Methodproject(a::Sequence, ::ParameterSpace, codomain_dest::VectorSpace, ::Type{T}=eltype(a))
Represent a
as a LinearOperator
from ParameterSpace
to codomain_dest
.
See also: project!
.
RadiiPolynomial.project
— Methodproject(a::Sequence, space_dest::VectorSpace, ::Type{T}=eltype(a))
Represent a
as a Sequence
in space_dest
.
See also: project!
.
RadiiPolynomial.project
— Methodproject(J::UniformScaling, domain_dest::VectorSpace, codomain_dest::VectorSpace, ::Type{T}=eltype(J))
Represent J
as a LinearOperator
from domain_dest
to codomain_dest
.
See also: project!
.
RadiiPolynomial.Multiplication
— TypeMultiplication{T<:Sequence{<:SequenceSpace}} <: SpecialOperator
Multiplication operator associated with a given Sequence
.
Field:
sequence :: T
Constructor:
Multiplication(::Sequence{<:SequenceSpace})
See also: *(::Sequence{<:SequenceSpace}, ::Sequence{<:SequenceSpace})
, mul!(::Sequence{<:SequenceSpace}, ::Sequence{<:SequenceSpace}, ::Sequence{<:SequenceSpace}, ::Number, ::Number)
, ^(::Sequence{<:SequenceSpace}, ::Int)
, project(::Multiplication, ::SequenceSpace, ::SequenceSpace)
, project!(::LinearOperator{<:SequenceSpace,<:SequenceSpace}, ::Multiplication)
and Multiplication
.
RadiiPolynomial.Multiplication
— Method(ℳ::Multiplication)(a::Sequence)
Compute the discrete convolution (associated with space(sequence(ℳ))
and space(a)
) of sequence(ℳ)
and a
; equivalent to sequence(ℳ) * a
.
See also: *(::Multiplication, ::Sequence)
, Multiplication
, *(::Sequence{<:SequenceSpace}, ::Sequence{<:SequenceSpace})
, mul!(::Sequence{<:SequenceSpace}, ::Sequence{<:SequenceSpace}, ::Sequence{<:SequenceSpace}, ::Number, ::Number)
and ^(::Sequence{<:SequenceSpace}, ::Int)
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator{<:SequenceSpace,<:SequenceSpace}, ℳ::Multiplication)
Represent ℳ
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project(::Multiplication, ::SequenceSpace, ::SequenceSpace)
and Multiplication
.
RadiiPolynomial.project
— Methodproject(ℳ::Multiplication, domain::SequenceSpace, codomain::SequenceSpace, ::Type{T}=eltype(sequence(ℳ)))
Represent ℳ
as a LinearOperator
from domain
to codomain
.
See also: project!(::LinearOperator{<:SequenceSpace,<:SequenceSpace}, ::Multiplication)
and Multiplication
.
RadiiPolynomial.Derivative
— TypeDerivative{T<:Union{Int,Tuple{Vararg{Int}}}} <: SpecialOperator
Generic derivative operator.
Field:
order :: T
Constructors:
Derivative(::Int)
Derivative(::Tuple{Vararg{Int}})
Derivative(order::Int...)
: equivalent toDerivative(order)
See also: differentiate
, differentiate!
, project(::Derivative, ::VectorSpace, ::VectorSpace)
and project!(::LinearOperator, ::Derivative)
.
Examples
julia> Derivative(1)
Derivative{Int64}(1)
julia> Derivative(1, 2)
Derivative{Tuple{Int64, Int64}}((1, 2))
RadiiPolynomial.Derivative
— Method(𝒟::Derivative)(a::AbstractSequence)
Compute the order(𝒟)
-th derivative of a
; equivalent to differentiate(a, order(𝒟))
.
See also: *(::Derivative, ::AbstractSequence)
, Derivative
, differentiate
and differentiate!
.
RadiiPolynomial.Integral
— TypeIntegral{T<:Union{Int,Tuple{Vararg{Int}}}} <: SpecialOperator
Generic integral operator.
Field:
order :: T
Constructors:
Integral(::Int)
Integral(::Tuple{Vararg{Int}})
Integral(order::Int...)
: equivalent toIntegral(order)
See also: integrate
, integrate!
, project(::Integral, ::VectorSpace, ::VectorSpace)
and project!(::LinearOperator, ::Integral)
.
Examples
julia> Integral(1)
Integral{Int64}(1)
julia> Integral(1, 2)
Integral{Tuple{Int64, Int64}}((1, 2))
RadiiPolynomial.Integral
— Method(ℐ::Integral)(a::AbstractSequence)
Compute the order(ℐ)
-th integral of a
; equivalent to integrate(a, order(ℐ))
.
See also: *(::Integral, ::AbstractSequence)
, Integral
, integrate
and integrate!
.
RadiiPolynomial.Laplacian
— TypeLaplacian <: SpecialOperator
Laplacian operator.
RadiiPolynomial.Laplacian
— Method(Δ::Laplacian)(a::AbstractSequence)
Compute the laplacian a
; equivalent to laplacian(a)
.
See also: *(::Laplacian, ::AbstractSequence)
, Laplacian
, laplacian
and laplacian!
.
RadiiPolynomial.differentiate
— Functiondifferentiate(a::Sequence, α=1)
Compute the α
-th derivative of a
.
See also: differentiate!
, Derivative
, *(::Derivative, ::Sequence)
and (::Derivative)(::Sequence)
.
RadiiPolynomial.differentiate!
— Functiondifferentiate!(c::Sequence, a::Sequence, α=1)
Compute the α
-th derivative of a
. The result is stored in c
by overwriting it.
See also: differentiate
, Derivative
, *(::Derivative, ::Sequence)
and (::Derivative)(::Sequence)
.
RadiiPolynomial.integrate
— Functionintegrate(a::Sequence, α=1)
Compute the α
-th integral of a
.
See also: integrate!
, Integral
, *(::Integral, ::Sequence)
and (::Integral)(::Sequence)
.
RadiiPolynomial.integrate!
— Functionintegrate!(c::Sequence, a::Sequence, α=1)
Compute the α
-th integral of a
. The result is stored in c
by overwriting it.
See also: integrate
, Integral
, *(::Integral, ::Sequence)
and (::Integral)(::Sequence)
.
RadiiPolynomial.laplacian!
— Methodlaplacian!(c::Sequence, a::Sequence)
Compute the laplacian a
. The result is stored in c
by overwriting it.
See also: laplacian
, Laplacian
, *(::Laplacian, ::Sequence)
and (::Laplacian)(::Sequence)
.
RadiiPolynomial.laplacian
— Methodlaplacian(a::Sequence)
Compute the laplacian of a
.
See also: laplacian!
, Laplacian
, *(::Laplacian, ::Sequence)
and (::Laplacian)(::Sequence)
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, 𝒟::Derivative)
Represent 𝒟
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project(::Derivative, ::VectorSpace, ::VectorSpace)
and Derivative
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, ℐ::Integral)
Represent ℐ
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project(::Integral, ::VectorSpace, ::VectorSpace)
and Integral
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, Δ::Laplacian)
Represent Δ
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project(::Laplacian, ::VectorSpace, ::VectorSpace)
and Laplacian
RadiiPolynomial.project
— Methodproject(𝒟::Derivative, domain::VectorSpace, codomain::VectorSpace, ::Type{T}=_coeftype(𝒟, domain, Float64))
Represent 𝒟
as a LinearOperator
from domain
to codomain
.
See also: project!(::LinearOperator, ::Derivative)
and Derivative
.
RadiiPolynomial.project
— Methodproject(ℐ::Integral, domain::VectorSpace, codomain::VectorSpace, ::Type{T}=_coeftype(ℐ, domain, Float64))
Represent ℐ
as a LinearOperator
from domain
to codomain
.
See also: project!(::LinearOperator, ::Integral)
and Integral
.
RadiiPolynomial.project
— Methodproject(Δ::Laplacian, domain::VectorSpace, codomain::VectorSpace, ::Type{T}=_coeftype(Δ, domain, Float64))
Represent Δ
as a LinearOperator
from domain
to codomain
.
See also: project!(::LinearOperator, ::Laplacian)
and Laplacian
.
RadiiPolynomial.Evaluation
— TypeEvaluation{T<:Union{Nothing,Number,Tuple{Vararg{Union{Nothing,Number}}}}} <: SpecialOperator
Generic evaluation operator. A value of nothing
indicates that no evaluation should be performed.
Field:
value :: T
Constructors:
Evaluation(::Union{Nothing,Number})
Evaluation(::Tuple{Vararg{Union{Nothing,Number}}})
Evaluation(value::Union{Number,Nothing}...)
: equivalent toEvaluation(value)
See also: evaluate
, evaluate!
, project(::Evaluation, ::VectorSpace, ::VectorSpace)
and project!(::LinearOperator, ::Evaluation)
.
Examples
julia> Evaluation(1.0)
Evaluation{Float64}(1.0)
julia> Evaluation(1.0, nothing, 2.0)
Evaluation{Tuple{Float64, Nothing, Float64}}((1.0, nothing, 2.0))
RadiiPolynomial.Evaluation
— Method(ℰ::Evaluation)(a::AbstractSequence)
Evaluate a
at value(ℰ)
; equivalent to evaluate(a, value(ℰ))
.
See also: *(::Evaluation, ::AbstractSequence)
, Evaluation
, (::AbstractSequence)(::Any, ::Vararg)
, evaluate
and evaluate!
.
RadiiPolynomial.evaluate!
— Methodevaluate!(c::Union{AbstractVector,Sequence}, a::Sequence, x)
Evaluate a
at x
. The result is stored in c
by overwriting it.
See also: (::Sequence)(::Any, ::Vararg)
, evaluate
, Evaluation
, *(::Evaluation, ::Sequence)
and (::Evaluation)(::Sequence)
.
RadiiPolynomial.evaluate
— Methodevaluate(a::Sequence, x)
Evaluate a
at x
.
See also: (::Sequence)(::Any, ::Vararg)
, evaluate!
, Evaluation
, *(::Evaluation, ::Sequence)
and (::Evaluation)(::Sequence)
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, ℰ::Evaluation)
Represent ℰ
as a LinearOperator
from domain
to codomain
. The result is stored in C
by overwriting it.
See also: project(::Evaluation, ::VectorSpace, ::VectorSpace)
and Evaluation
.
RadiiPolynomial.project
— Methodproject(ℰ::Evaluation, domain::VectorSpace, codomain::VectorSpace, ::Type{T}=_coeftype(ℰ, domain, Float64))
Represent ℰ
as a LinearOperator
from domain
to codomain
.
See also: project!(::LinearOperator, ::Evaluation)
and Evaluation
.
RadiiPolynomial.Scale
— TypeScale{T<:Union{Number,Tuple{Vararg{Number}}}} <: SpecialOperator
Generic scale operator.
Field:
value :: T
Constructors:
Scale(::Number)
Scale(::Tuple{Vararg{Number}})
Scale(value::Number...)
: equivalent toScale(value)
See also: scale
, scale!
, project(::Scale, ::VectorSpace, ::VectorSpace)
and project!(::LinearOperator, ::Scale)
.
Examples
julia> Scale(1.0)
Scale{Float64}(1.0)
julia> Scale(1.0, 2.0)
Scale{Tuple{Float64, Float64}}((1.0, 2.0))
RadiiPolynomial.Scale
— Method(𝒮::Scale)(a::AbstractSequence)
Scale a
by a factor value(𝒮)
; equivalent to scale(a, value(𝒮))
.
See also: *(::Scale, ::AbstractSequence)
, Scale
, scale
and scale!
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, 𝒮::Scale)
Represent 𝒮
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project(::Scale, ::VectorSpace, ::VectorSpace)
and Scale
RadiiPolynomial.project
— Methodproject(𝒮::Scale, domain::VectorSpace, codomain::VectorSpace, ::Type{T}=_coeftype(𝒮, domain, Float64))
Represent 𝒮
as a LinearOperator
from domain
to codomain
.
See also: project!(::LinearOperator, ::Scale)
and Scale
RadiiPolynomial.scale!
— Methodscale!(c::Sequence, a::Sequence, γ)
Scale a
by a factor γ
. The result is stored in c
by overwriting it.
See also: scale
, Scale
, *(::Scale, ::Sequence)
and (::Scale)(::Sequence)
.
RadiiPolynomial.scale
— Methodscale(a::Sequence, γ)
Scale a
by a factor γ
.
See also: scale!
, Scale
, *(::Scale, ::Sequence)
and (::Scale)(::Sequence)
.
RadiiPolynomial.Shift
— TypeShift{T<:Union{Number,Tuple{Vararg{Number}}}} <: SpecialOperator
Generic shift operator.
Field:
value :: T
Constructors:
Shift(::Number)
Shift(::Tuple{Vararg{Number}})
Shift(value::Number...)
: equivalent toShift(value)
See also: shift
, shift!
, project(::Shift, ::VectorSpace, ::VectorSpace)
and project!(::LinearOperator, ::Shift)
.
Examples
julia> Shift(1.0)
Shift{Float64}(1.0)
julia> Shift(1.0, 2.0)
Shift{Tuple{Float64, Float64}}((1.0, 2.0))
RadiiPolynomial.Shift
— Method(𝒮::Shift)(a::AbstractSequence)
Shift a
by value(𝒮)
; equivalent to shift(a, value(𝒮))
.
See also: *(::Shift, ::AbstractSequence)
, Shift
, shift
and shift!
.
RadiiPolynomial.project!
— Methodproject!(C::LinearOperator, 𝒮::Shift)
Represent 𝒮
as a LinearOperator
from domain(C)
to codomain(C)
. The result is stored in C
by overwriting it.
See also: project(::Shift, ::VectorSpace, ::VectorSpace)
and Shift
.
RadiiPolynomial.project
— Methodproject(𝒮::Shift, domain::VectorSpace, codomain::VectorSpace, ::Type{T}=_coeftype(𝒮, domain, Float64))
Represent 𝒮
as a LinearOperator
from domain
to codomain
.
See also: project!(::LinearOperator, ::Shift)
and Shift
.
RadiiPolynomial.shift!
— Methodshift!(c::Sequence, a::Sequence, τ)
Shift a
by τ
. The result is stored in c
by overwriting it.
See also: shift
, Shift
, *(::Shift, ::Sequence)
and (::Shift)(::Sequence)
.
RadiiPolynomial.shift
— Methodshift(a::Sequence, τ)
Shift a
by τ
.
See also: shift!
, Shift
, *(::Shift, ::Sequence)
and (::Shift)(::Sequence)
.