Reference/API
All polynomials have the following functionality. In some cases, there is not a direct function call and therefore the polynomials have to be converted to the standard Polynomial
type before continuing.
Base.chop
Base.gcd
Base.length
Base.one
Base.size
Base.truncate
Base.zero
Polynomials.chop!
Polynomials.coeffs
Polynomials.companion
Polynomials.degree
Polynomials.derivative
Polynomials.domain
Polynomials.fit
Polynomials.fromroots
Polynomials.integrate
Polynomials.mapdomain
Polynomials.roots
Polynomials.truncate!
Polynomials.vander
Polynomials.variable
Inspection
Polynomials.coeffs
— Functioncoeffs(::AbstractPolynomial)
Return the coefficient vector [a_0, a_1, ..., a_n]
of a polynomial.
Polynomials.degree
— Functiondegree(::AbstractPolynomial)
Return the degree of the polynomial, i.e. the highest exponent in the polynomial that has a nonzero coefficient. The degree of the zero polynomial is defined to be -1.
Base.length
— Functionlength(::AbstractPolynomial)
The length of the polynomial.
Base.size
— Functionsize(::AbstractPolynomial, [i])
Returns the size of the polynomials coefficients, along axis i
if provided.
Polynomials.domain
— Functiondomain(::Type{<:AbstractPolynomial})
Returns the domain of the polynomial.
Polynomials.mapdomain
— Functionmapdomain(::Type{<:AbstractPolynomial}, x::AbstractArray)
mapdomain(::AbstractPolynomial, x::AbstractArray)
Given values of x that are assumed to be unbounded (-∞, ∞), return values rescaled to the domain of the given polynomial.
Examples
julia> using Polynomials
julia> x = -10:10
-10:10
julia> extrema(mapdomain(ChebyshevT, x))
(-1.0, 1.0)
Base.chop
— Functionchop(::AbstractPolynomial{T};
rtol::Real = Base.rtoldefault(real(T)), atol::Real = 0))
Removes any leading coefficients that are approximately 0 (using rtol
and atol
). Returns a polynomial whose degree will guaranteed to be equal to or less than the given polynomial's.
Polynomials.chop!
— Functionchop!(::AbstractPolynomial{T};
rtol::Real = Base.rtoldefault(real(T)), atol::Real = 0))
In-place version of chop
Base.truncate
— Functiontruncate(::AbstractPolynomial{T};
rtol::Real = Base.rtoldefault(real(T)), atol::Real = 0)
Rounds off coefficients close to zero, as determined by rtol
and atol
, and then chops any leading zeros. Returns a new polynomial.
Polynomials.truncate!
— Functiontruncate!(::AbstractPolynomial{T};
rtol::Real = Base.rtoldefault(real(T)), atol::Real = 0)
In-place version of truncate
Arithmetic
All AbstractPolynomials
have basic arithmetic operations defined on them (+
, -
, *
, /
, ÷
, %
, ==
).
julia> p = Polynomial([1, 2])
Polynomial(1 + 2*x)
julia> q = Polynomial([1, 0, -1])
Polynomial(1 - x^2)
julia> 2p
Polynomial(2 + 4*x)
julia> 2 + p
Polynomial(3 + 2*x)
julia> p - q
Polynomial(2*x + x^2)
julia> p * q
Polynomial(1 + 2*x - x^2 - 2*x^3)
julia> q / 2
Polynomial(0.5 - 0.5*x^2)
Base.gcd
— Functiongcd(a::AbstractPolynomial, b::AbstractPolynomial; atol::Real=0, rtol::Real=Base.rtoldefault)
Find the greatest common denominator of two polynomials recursively using Euclid's algorithm.
Examples
julia> using Polynomials
julia> gcd(fromroots([1, 1, 2]), fromroots([1, 2, 3]))
Polynomial(4.0 - 6.0*x + 2.0*x^2)
gcd(p1::StandardBasisPolynomial, p2::StandardBasisPolynomial; method=:eculidean, kwargs...)
Find the greatest common divisor.
By default, uses the Euclidean division algorithm (method=:euclidean
), which is susceptible to floating point issues.
Passing method=:noda_sasaki
uses scaling to circumvent some of these.
Passing method=:numerical
will call the internal method NGCD.ngcd
for the numerical gcd. See the help page of Polynomials.NGCD.ngcd
for details.
Mathematical Functions
Base.zero
— Functionzero(::Type{<:AbstractPolynomial})
zero(::AbstractPolynomial)
Returns a representation of 0 as the given polynomial.
Base.one
— Functionone(::Type{<:AbstractPolynomial})
one(::AbstractPolynomial)
Returns a representation of 1 as the given polynomial.
Polynomials.variable
— Functionvariable(var=:x)
variable(::Type{<:AbstractPolynomial}, var=:x)
variable(p::AbstractPolynomial, var=p.var)
Return the monomial x
in the indicated polynomial basis. If no type is give, will default to Polynomial
. Equivalent to P(var)
.
Examples
julia> using Polynomials
julia> x = variable()
Polynomial(x)
julia> p = 100 + 24x - 3x^2
Polynomial(100 + 24*x - 3*x^2)
julia> roots((x - 3) * (x + 2))
2-element Array{Float64,1}:
-2.0
3.0
Polynomials.fromroots
— Functionfromroots(::AbstractVector{<:Number}; var=:x)
fromroots(::Type{<:AbstractPolynomial}, ::AbstractVector{<:Number}; var=:x)
Construct a polynomial of the given type given the roots. If no type is given, defaults to Polynomial
.
Examples
julia> using Polynomials
julia> r = [3, 2]; # (x - 3)(x - 2)
julia> fromroots(r)
Polynomial(6 - 5*x + x^2)
fromroots(::AbstractMatrix{<:Number}; var=:x)
fromroots(::Type{<:AbstractPolynomial}, ::AbstractMatrix{<:Number}; var=:x)
Construct a polynomial of the given type using the eigenvalues of the given matrix as the roots. If no type is given, defaults to Polynomial
.
Examples
julia> using Polynomials
julia> A = [1 2; 3 4]; # (x - 5.37228)(x + 0.37228)
julia> fromroots(A)
Polynomial(-1.9999999999999998 - 5.0*x + 1.0*x^2)
Polynomials.roots
— Functionroots(::AbstractPolynomial; kwargs...)
Returns the roots of the given polynomial. This is calculated via the eigenvalues of the companion matrix. The kwargs
are passed to the LinearAlgeebra.eigvals
call.
The [PolynomialRoots.jl](https://github.com/giordano/PolynomialRoots.jl) package provides an alternative that is a bit faster and a bit more accurate; the [FastPolynomialRoots](https://github.com/andreasnoack/FastPolynomialRoots.jl) provides an interface to FORTRAN code implementing an algorithm that can handle very large polynomials (it is `O(n^2)` not `O(n^3)`. the [AMRVW.jl](https://github.com/jverzani/AMRVW.jl) package implements the algorithm in Julia, allowing the use of other number types.
roots(p)
Compute the roots of the Laurent polynomial p
.
The roots of a function (Laurent polynomial in this case) a(z)
are the values of z
for which the function vanishes. A Laurent polynomial $a(z) = a_m z^m + a_{m+1} z^{m+1} + ... + a_{-1} z^{-1} + a_0 + a_1 z + ... + a_{n-1} z^{n-1} + a_n z^n$ can equivalently be viewed as a rational function with a multiple singularity (pole) at the origin. The roots are then the roots of the numerator polynomial. For example, $a(z) = 1/z + 2 + z$ can be written as $a(z) = (1+2z+z^2) / z$ and the roots of a
are the roots of $1+2z+z^2$.
Example
julia> p = LaurentPolynomial([24,10,-15,0,1],-2:1,:z)
LaurentPolynomial(24*z⁻² + 10*z⁻¹ - 15 + z²)
julia> roots(a)
4-element Array{Float64,1}:
-3.999999999999999
-0.9999999999999994
1.9999999999999998
2.9999999999999982
Polynomials.derivative
— Functionderivative(::AbstractPolynomial, order::Int = 1)
Returns a polynomial that is the order
th derivative of the given polynomial. order
must be non-negative.
Polynomials.integrate
— Functionintegrate(::AbstractPolynomial, C=0)
Returns the indefinite integral of the polynomial with constant C
.
integrate(::AbstractPolynomial, a, b)
Compute the definite integral of the given polynomial from a
to b
. Will throw an error if either a
or b
are out of the polynomial's domain.
Polynomials.fit
— Functionfit(x, y, deg=length(x) - 1; [weights], var=:x)
fit(::Type{<:AbstractPolynomial}, x, y, deg=length(x)-1; [weights], var=:x)
Fit the given data as a polynomial type with the given degree. Uses linear least squares. When weights are given, as either a Number
, Vector
or Matrix
, will use weighted linear least squares. The default polynomial type is Polynomial
. This will automatically scale your data to the domain
of the polynomial type using mapdomain
Polynomials.companion
— Functioncompanion(::AbstractPolynomial)
Return the companion matrix for the given polynomial.
References
Polynomials.vander
— Functionvander(::Type{AbstractPolynomial}, x::AbstractVector, deg::Integer)
Calculate the psuedo-Vandermonde matrix of the given polynomial type with the given degree.
References
Plotting
Polynomials can be plotted directly using Plots.jl.
plot(::AbstractPolynomial; kwds...)
will automatically determine a range based on the critical points (roots, extrema and points of inflection).
plot(::AbstractPolynomial, a, b; kwds...)
will plot the polynomial within the range [a, b]
.
Example: The Polynomials.jl logo
using Plots, Polynomials
# T1, T2, T3, and T4:
chebs = [
ChebyshevT([0, 1]),
ChebyshevT([0, 0, 1]),
ChebyshevT([0, 0, 0, 1]),
ChebyshevT([0, 0, 0, 0, 1]),
]
colors = ["#4063D8", "#389826", "#CB3C33", "#9558B2"]
itr = zip(chebs, colors)
(cheb,col), state = iterate(itr)
p = plot(cheb, c=col, lw=5, legend=false, label="")
for (cheb, col) in Base.Iterators.rest(itr, state)
plot!(cheb, c=col, lw=5)
end