Library
TaylorSeries
— Module.TaylorSeries
A Julia package for Taylor expansions in one or more independent variables.
The basic constructors are Taylor1
and TaylorN
; see also HomogeneousPolynomial
.
Types
TaylorSeries.Taylor1
— Type.Taylor1{T<:Number} <: AbstractSeries{T}
DataType for polynomial expansions in one independent variable.
Fields:
coeffs :: Array{T,1}
Expansion coefficients; the $i$-th component is the coefficient of degree $i-1$ of the expansion.order :: Int64
Maximum order (degree) of the polynomial.
Note that Taylor1
variables are callable. For more information, see evaluate
.
HomogeneousPolynomial{T<:Number} <: AbstractSeries{T}
DataType for homogenous polynomials in many (>1) independent variables.
Fields:
coeffs :: Array{T,1}
Expansion coefficients of the homogeneous
polynomial; the $i$-th component is related to a monomial, where the degrees of the independent variables are specified by coeff_table[order+1][i]
.
order :: Int
order (degree) of the homogenous polynomial.
Note that HomogeneousPolynomial
variables are callable. For more information, see evaluate
.
TaylorSeries.TaylorN
— Type.TaylorN{T<:Number} <: AbstractSeries{T}
DataType for polynomial expansions in many (>1) independent variables.
Fields:
coeffs :: Array{HomogeneousPolynomial{T},1}
Vector containing the
HomogeneousPolynomial
entries. The $i$-th component corresponds to the homogeneous polynomial of degree $i-1$.
order :: Int
maximum order of the polynomial expansion.
Note that TaylorN
variables are callable. For more information, see evaluate
.
TaylorSeries.AbstractSeries
— Type.AbstractSeries{T<:Number} <: Number
Parameterized abstract type for Taylor1
, HomogeneousPolynomial
and TaylorN
.
TaylorSeries.ParamsTaylorN
— Type.ParamsTaylorN
DataType holding the current parameters for TaylorN
and HomogeneousPolynomial
.
Fields:
order :: Int
Order (degree) of the polynomialsnum_vars :: Int
Number of variablesvariable_names :: Array{String,1}
Name of the variables
These parameters can be changed using set_params_TaylorN(order, numVars)
.
Functions and methods
TaylorSeries.Taylor1
— Method.Taylor1([T::Type=Float64], [order::Int=1])
Shortcut to define the independent variable of a Taylor1{T}
polynomial of given order
. The default type for T
is Float64
.
julia> Taylor1(16)
1.0 t + 𝒪(t¹⁷)
julia> Taylor1(Rational{Int}, 4)
1//1 t + 𝒪(t⁵)
TaylorSeries.HomogeneousPolynomial
— Method.HomogeneousPolynomial([T::Type=Float64], nv::Int])
Shortcut to define the nv
-th independent HomogeneousPolynomial{T}
. The default type for T
is Float64
.
julia> HomogeneousPolynomial(1)
1.0 x₁
julia> HomogeneousPolynomial(Rational{Int}, 2)
1//1 x₂
TaylorSeries.TaylorN
— Method.TaylorN([T::Type=Float64], nv::Int; [order::Int=get_order()])
Shortcut to define the nv
-th independent TaylorN{T}
variable as a polynomial. The order is defined through the keyword parameter order
, whose default corresponds to get_order()
. The default of type for T
is Float64
.
julia> TaylorN(1)
1.0 x₁ + 𝒪(‖x‖⁷)
julia> TaylorN(Rational{Int},2)
1//1 x₂ + 𝒪(‖x‖⁷)
TaylorSeries.set_variables
— Function.set_variables([T::Type], names::String; [order=get_order(), numvars=-1])
Return a TaylorN{T}
vector with each entry representing an independent variable. names
defines the output for each variable (separated by a space). The default type T
is Float64
, and the default for order
is the one defined globally. Changing the order
or numvars
resets the hash_tables.
If numvars
is not specified, it is inferred from names
. If only one variable name is defined and numvars>1
, it uses this name with subscripts for the different variables.
julia> set_variables(Int, "x y z", order=4)
3-element Array{TaylorSeries.TaylorN{Int64},1}:
1 x + 𝒪(‖x‖⁵)
1 y + 𝒪(‖x‖⁵)
1 z + 𝒪(‖x‖⁵)
julia> set_variables("α", numvars=2)
2-element Array{TaylorSeries.TaylorN{Float64},1}:
1.0 α₁ + 𝒪(‖x‖⁵)
1.0 α₂ + 𝒪(‖x‖⁵)
julia> set_variables("x", order=6, numvars=2)
2-element Array{TaylorSeries.TaylorN{Float64},1}:
1.0 x₁ + 𝒪(‖x‖⁷)
1.0 x₂ + 𝒪(‖x‖⁷)
TaylorSeries.get_variables
— Function.get_variables(;order=get_order())
Return a TaylorN
vector with each entry representing an independent variable. It takes the default _params_TaylorN_
values if set_variables
hasn't been changed with the exception that order
can be explicitely established by the user without changing internal values for num_vars
or variable_names
.
TaylorSeries.show_params_TaylorN
— Function.show_params_TaylorN()
Display the current parameters for TaylorN
and HomogeneousPolynomial
types.
TaylorSeries.show_monomials
— Function.show_monomials(ord::Int) --> nothing
List the indices and corresponding of a HomogeneousPolynomial
of degree ord
.
TaylorSeries.getcoeff
— Function.getcoeff(a, n)
Return the coefficient of order n::Int
of a a::Taylor1
polynomial.
getcoeff(a, v)
Return the coefficient of a::HomogeneousPolynomial
, specified by v::Array{Int,1}
which has the indices of the specific monomial.
getcoeff(a, v)
Return the coefficient of a::TaylorN
, specified by v::Array{Int,1}
which has the indices of the specific monomial.
TaylorSeries.evaluate
— Function.evaluate(a, [dx])
Evaluate a Taylor1
polynomial using Horner's rule (hand coded). If dx
is ommitted, its value is considered as zero. Note that a Taylor1
polynomial a
may also be evaluated by calling it as a function; that is, the syntax a(dx)
is equivalent to evaluate(a,dx)
, and a()
is equivalent to evaluate(a)
.
evaluate(x, δt)
Evaluates each element of x::Array{Taylor1{T},1}
, representing the dependent variables of an ODE, at time δt. Note that an array x
of Taylor1
polynomials may also be evaluated by calling it as a function; that is, the syntax x(δt)
is equivalent to evaluate(x, δt)
, and x()
is equivalent to evaluate(x)
.
evaluate(a, x)
Substitute x::Taylor1
as independent variable in a a::Taylor1
polynomial. Note that the syntax a(x)
is equivalent to evaluate(a, x)
.
evaluate(a, [vals])
Evaluate a HomogeneousPolynomial
polynomial at vals
. If vals
is ommitted, it's evaluated at zero. Note that the syntax a(vals)
is equivalent to evaluate(a, vals)
; and a()
is equivalent to evaluate(a)
.
evaluate(a, [vals])
Evaluate the TaylorN
polynomial a
at vals
. If vals
is ommitted, it's evaluated at zero. Note that the syntax a(vals)
is equivalent to evaluate(a, vals)
; and a()
is equivalent to evaluate(a)
.
TaylorSeries.evaluate!
— Function.evaluate!(x, δt, x0)
Evaluates each element of x::Array{Taylor1{T},1}
, representing the Taylor expansion for the dependent variables of an ODE at time δt
. It updates the vector x0
with the computed values.
TaylorSeries.taylor_expand
— Function.taylor_expand(f, x0; order)
Computes the Taylor expansion of the function f
around the point x0
.
If x0
is a scalar, a Taylor1
expansion will be returned. If x0
is a vector, a TaylorN
expansion will be computed. If the dimension of x0 (length(x0)
) is different from the variables set for TaylorN
(get_numvars()
), an AssertionError
will be thrown.
TaylorSeries.update!
— Function.update!(a, x0)
Takes a <: Union{Taylo1,TaylorN}
and expands it around the coordinate x0
.
TaylorSeries.derivative
— Function.derivative(a)
Return the Taylor1
polynomial of the differential of a::Taylor1
. The last coefficient is set to zero.
derivative(a, n)
Compute recursively the Taylor1
polynomial of the n-th derivative of a::Taylor1
.
derivative(n, a)
Return the value of the n
-th derivative of the polynomial a
.
derivative(a, r)
Partial differentiation of a::HomogeneousPolynomial
series with respect to the r
-th variable.
derivative(a, [r=1])
Partial differentiation of a::TaylorN
series with respect to the r
-th variable.
TaylorSeries.integrate
— Function.integrate(a, [x])
Return the integral of a::Taylor1
. The constant of integration (0-th order coefficient) is set to x
, which is zero if ommitted.
integrate(a, r)
Integrate the a::HomogeneousPolynomial
with respect to the r
-th variable. The returned HomogeneousPolynomial
has no added constant of integration. If the order of a corresponds to get_order()
, a zero HomogeneousPolynomial
of 0-th order is returned.
integrate(a, r, [x0])
Integrate the a::TaylorN
series with respect to the r
-th variable, where x0
the integration constant and must be independent of the r
-th variable; if x0
is ommitted, it is taken as zero.
Base.LinAlg.gradient
— Function. gradient(f)
∇(f)
Compute the gradient of the polynomial f::TaylorN
.
TaylorSeries.jacobian
— Function. jacobian(vf)
jacobian(vf, [vals])
Compute the jacobian matrix of vf
, a vector of TaylorN
polynomials, evaluated at the vector vals
. If vals
is ommited, it is evaluated at zero.
TaylorSeries.jacobian!
— Function. jacobian!(jac, vf)
jacobian!(jac, vf, [vals])
Compute the jacobian matrix of vf
, a vector of TaylorN
polynomials evaluated at the vector vals
, and write results to jac
. If vals
is ommited, it is evaluated at zero.
TaylorSeries.hessian
— Function. hessian(f)
hessian(f, [vals])
Return the hessian matrix (jacobian of the gradient) of f::TaylorN
, evaluated at the vector vals
. If vals
is ommited, it is evaluated at zero.
TaylorSeries.hessian!
— Function. hessian!(hes, f)
hessian!(hes, f, [vals])
Return the hessian matrix (jacobian of the gradient) of f::TaylorN
, evaluated at the vector vals
, and write results to hes
. If vals
is ommited, it is evaluated at zero.
TaylorSeries.inverse
— Function.inverse(f)
Return the Taylor expansion of $f^{-1}(t)$, of order N = f.order
, for f::Taylor1
polynomial if the first coefficient of f
is zero. Otherwise, an ArgumentError
is thrown.
The algorithm implements Lagrange inversion at $t=0$ if $f(0)=0$:
Base.abs
— Function.abs(a)
Returns a
if constant_term(a) > 0
and -a
if constant_term(a) < 0
for a <:Union{Taylor1,TaylorN}
. Notice that typeof(abs(a)) <: AbstractSeries
.
Base.LinAlg.norm
— Function.norm(x::AbstractSeries, p::Real)
Returns the p-norm of an x::AbstractSeries
, defined by
which returns a non-negative number.
Base.isapprox
— Function.isapprox(x::AbstractSeries, y::AbstractSeries;
rtol::Real=sqrt(eps), atol::Real=0, nans::Bool=false)
Inexact equality comparison between polynomials: returns true
if norm(x-y,1) <= atol + rtol*max(norm(x,1), norm(y,1))
, where x
and y
are polynomials. For more details, see Base.isapprox
.
Base.isfinite
— Function.isfinite(x::AbstractSeries) -> Bool
Test whether the coefficients of the polynomial x
are finite.
TaylorSeries.displayBigO
— Function.displayBigO(d::Bool) --> nothing
Set/unset displaying of the big 𝒪 notation in the output of Taylor1
and TaylorN
polynomials. The initial value is true
.
Internals
TaylorSeries.generate_tables
— Function.generate_tables(num_vars, order)
Return the hash tables coeff_table
, index_table
, size_table
and pos_table
. Internally, these are treated as const
.
Hash tables
coeff_table :: Array{Array{Array{Int64,1},1},1}
The $i+1$-th component contains a vector with the vectors of all the possible combinations of monomials of a HomogeneousPolynomial
of order $i$.
index_table :: Array{Array{Int64,1},1}
The $i+1$-th component contains a vector of (hashed) indices that represent the distinct monomials of a HomogeneousPolynomial
of order (degree) $i$.
size_table :: Array{Int64,1}
The $i+1$-th component contains the number of distinct monomials of the HomogeneousPolynomial
of order $i$, equivalent to length(coeff_table[i])
.
pos_table :: Array{Dict{Int64,Int64},1}
The $i+1$-th component maps the hash index to the (lexicographic) position of the corresponding monomial in coeffs_table
.
TaylorSeries.generate_index_vectors
— Function.generate_index_vectors(num_vars, degree)
Return a vector of index vectors with num_vars
(number of variables) and degree.
TaylorSeries.in_base
— Function.in_base(order, v)
Convert vector v
of non-negative integers to base order+1
.
TaylorSeries.make_inverse_dict
— Function.make_inverse_dict(v)
Return a Dict with the enumeration of v
: the elements of v
point to the corresponding index.
It is used to construct pos_table
from index_table
.
TaylorSeries.resize_coeffs1!
— Function.resize_coeffs1!{T<Number}(coeffs::Array{T,1}, order::Int)
If the length of coeffs
is smaller than order+1
, it resizes coeffs
appropriately filling it with zeros.
TaylorSeries.resize_coeffsHP!
— Function.resize_coeffsHP!{T<Number}(coeffs::Array{T,1}, order::Int)
If the length of coeffs
is smaller than the number of coefficients correspondinf to order
(given by size_table[order+1]
), it resizes coeffs
appropriately filling it with zeros.
TaylorSeries.constant_term
— Function.constant_term(a)
Return the constant value (zero order coefficient) for Taylor1
and TaylorN
.
TaylorSeries.mul!
— Function.mul!(c, a, b, k::Int) --> nothing
Update the k
-th expansion coefficient c[k]
of c = a * b
, where all c
, a
, and b
are either Taylor1
or TaylorN
.
The coefficients are given by
TaylorSeries.mul!
— Method.mul!(c, a, b) --> nothing
Return c = a*b
with no allocation; all arguments are HomogeneousPolynomial
.
TaylorSeries.div!
— Function.div!(c, a, b, k::Int, ordfact::Int=0)
Compute the k-th
expansion coefficient c[k]
of c = a / b
, where all c
, a
and b
are either Taylor1
or TaylorN
.
The coefficients are given by
For Taylor1
polynomials, ordfact
is the order of the factorized term of the denominator.
TaylorSeries.pow!
— Function.pow!(c, a, r::Real, k::Int, k0::Int=0)
Update the k-th
expansion coefficient c[k]
of c = a^r
, for both c
and a
either Taylor1
or TaylorN
.
The coefficients are given by
For Taylor1
polynomials, k0
is the order of the first non-zero coefficient of a
.
TaylorSeries.square
— Function.square(a::AbstractSeries) --> typeof(a)
Return a^2
; see TaylorSeries.sqr!
.
TaylorSeries.sqr!
— Function.sqr!(c, a, k::Int) --> nothing
Update the k-th
expansion coefficient c[k]
of c = a^2
, for both c
and a
either Taylor1
or TaylorN
.
The coefficients are given by
TaylorSeries.sqr!
— Method.sqr!(c, a)
Return c = a*a
with no allocation; all parameters are HomogeneousPolynomial
.
TaylorSeries.sqrt!
— Function.sqrt!(c, a, k::Int, k0::Int=0)
Compute the k-th
expansion coefficient c[k]
of c = sqrt(a)
for bothc
and a
either Taylor1
or TaylorN
.
The coefficients are given by
For Taylor1
polynomials, k0
is the order of the first non-zero coefficient, which must be even.
TaylorSeries.exp!
— Function.exp!(c, a, k) --> nothing
Update the k-th
expansion coefficient c[k+1]
of c = exp(a)
for both c
and a
either Taylor1
or TaylorN
.
The coefficients are given by
TaylorSeries.log!
— Function.log!(c, a, k) --> nothing
Update the k-th
expansion coefficient c[k+1]
of c = log(a)
for both c
and a
either Taylor1
or TaylorN
.
The coefficients are given by
TaylorSeries.sincos!
— Function.sincos!(s, c, a, k) --> nothing
Update the k-th
expansion coefficients s[k+1]
and c[k+1]
of s = sin(a)
and c = cos(a)
simultaneously, for s
, c
and a
either Taylor1
or TaylorN
.
The coefficients are given by
TaylorSeries.tan!
— Function.tan!(c, a, p, k::Int) --> nothing
Update the k-th
expansion coefficients c[k+1]
of c = tan(a)
, for c
and a
either Taylor1
or TaylorN
; p = c^2
and is passed as an argument for efficiency.
The coefficients are given by
TaylorSeries.asin!
— Function.asin!(c, a, r, k)
Update the k-th
expansion coefficients c[k+1]
of c = asin(a)
, for c
and a
either Taylor1
or TaylorN
; r = sqrt(1-c^2)
and is passed as an argument for efficiency.
TaylorSeries.acos!
— Function.acos!(c, a, r, k)
Update the k-th
expansion coefficients c[k+1]
of c = acos(a)
, for c
and a
either Taylor1
or TaylorN
; r = sqrt(1-c^2)
and is passed as an argument for efficiency.
TaylorSeries.atan!
— Function.atan!(c, a, r, k)
Update the k-th
expansion coefficients c[k+1]
of c = atan(a)
, for c
and a
either Taylor1
or TaylorN
; r = 1+a^2
and is passed as an argument for efficiency.
TaylorSeries.sinhcosh!
— Function.sinhcosh!(s, c, a, k)
Update the k-th
expansion coefficients s[k+1]
and c[k+1]
of s = sinh(a)
and c = cosh(a)
simultaneously, for s
, c
and a
either Taylor1
or TaylorN
.
The coefficients are given by
TaylorSeries.tanh!
— Function.tanh!(c, a, p, k)
Update the k-th
expansion coefficients c[k+1]
of c = tanh(a)
, for c
and a
either Taylor1
or TaylorN
; p = a^2
and is passed as an argument for efficiency.
Base.LinAlg.A_mul_B!
— Function.A_mul_B!(Y, A, B)
Multiply A*B and save the result in Y.
TaylorSeries.derivative!
— Function.derivative!(res, a) --> nothing
In-place version of derivative
. Compute the Taylor1
polynomial of the differential of a::Taylor1
and save it into res
. The last coefficient is set to zero.
derivative!(p, a, k) --> nothing
Update in-place the k-th
expansion coefficient p[k]
of p = derivative(a)
for both p
and a
Taylor1
.
The coefficients are given by
Index
TaylorSeries.AbstractSeries
TaylorSeries.HomogeneousPolynomial
TaylorSeries.HomogeneousPolynomial
TaylorSeries.ParamsTaylorN
TaylorSeries.Taylor1
TaylorSeries.Taylor1
TaylorSeries.TaylorN
TaylorSeries.TaylorN
Base.LinAlg.A_mul_B!
Base.LinAlg.gradient
Base.LinAlg.norm
Base.abs
Base.isapprox
Base.isfinite
TaylorSeries.acos!
TaylorSeries.asin!
TaylorSeries.atan!
TaylorSeries.constant_term
TaylorSeries.derivative
TaylorSeries.derivative!
TaylorSeries.displayBigO
TaylorSeries.div!
TaylorSeries.evaluate
TaylorSeries.evaluate!
TaylorSeries.exp!
TaylorSeries.generate_index_vectors
TaylorSeries.generate_tables
TaylorSeries.get_variables
TaylorSeries.getcoeff
TaylorSeries.hessian
TaylorSeries.hessian!
TaylorSeries.in_base
TaylorSeries.integrate
TaylorSeries.inverse
TaylorSeries.jacobian
TaylorSeries.jacobian!
TaylorSeries.log!
TaylorSeries.make_inverse_dict
TaylorSeries.mul!
TaylorSeries.mul!
TaylorSeries.pow!
TaylorSeries.resize_coeffs1!
TaylorSeries.resize_coeffsHP!
TaylorSeries.set_variables
TaylorSeries.show_monomials
TaylorSeries.show_params_TaylorN
TaylorSeries.sincos!
TaylorSeries.sinhcosh!
TaylorSeries.sqr!
TaylorSeries.sqr!
TaylorSeries.sqrt!
TaylorSeries.square
TaylorSeries.tan!
TaylorSeries.tanh!
TaylorSeries.taylor_expand
TaylorSeries.update!