API

Library


TaylorSeriesModule.
TaylorSeries

A Julia package for Taylor expansions in one or more independent variables.

The basic constructors are Taylor1 and TaylorN; see also HomogeneousPolynomial.

source

Types

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.

source
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.

source
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.

source
AbstractSeries{T<:Number} <: Number

Parameterized abstract type for Taylor1, HomogeneousPolynomial and TaylorN.

source
ParamsTaylorN

DataType holding the current parameters for TaylorN and HomogeneousPolynomial.

Fields:

  • order :: Int Order (degree) of the polynomials

  • num_vars :: Int Number of variables

  • variable_names :: Array{String,1} Name of the variables

These parameters can be changed using set_params_TaylorN(order, numVars).

source

Functions and methods

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⁵)
source
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₂
source
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‖⁷)
source
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‖⁷)
source
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.

source
show_params_TaylorN()

Display the current parameters for TaylorN and HomogeneousPolynomial types.

source
show_monomials(ord::Int) --> nothing

List the indices and corresponding of a HomogeneousPolynomial of degree ord.

source
TaylorSeries.getcoeffFunction.
getcoeff(a, n)

Return the coefficient of order n::Int of a a::Taylor1 polynomial.

source
getcoeff(a, v)

Return the coefficient of a::HomogeneousPolynomial, specified by v::Array{Int,1} which has the indices of the specific monomial.

source
getcoeff(a, v)

Return the coefficient of a::TaylorN, specified by v::Array{Int,1} which has the indices of the specific monomial.

source
TaylorSeries.evaluateFunction.
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).

source
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).

source
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).

source
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).

source
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).

source
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.

source
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.

source
TaylorSeries.update!Function.
update!(a, x0)

Takes a <: Union{Taylo1,TaylorN} and expands it around the coordinate x0.

source
derivative(a)

Return the Taylor1 polynomial of the differential of a::Taylor1. The last coefficient is set to zero.

source
derivative(a, n)

Compute recursively the Taylor1 polynomial of the n-th derivative of a::Taylor1.

source
derivative(n, a)

Return the value of the n-th derivative of the polynomial a.

source
derivative(a, r)

Partial differentiation of a::HomogeneousPolynomial series with respect to the r-th variable.

source
derivative(a, [r=1])

Partial differentiation of a::TaylorN series with respect to the r-th variable.

source
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.

source
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.

source
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.

source
Base.LinAlg.gradientFunction.
    gradient(f)
    ∇(f)

Compute the gradient of the polynomial f::TaylorN.

source
TaylorSeries.jacobianFunction.
    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.

source
    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.

source
TaylorSeries.hessianFunction.
    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.

source
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.

source
TaylorSeries.inverseFunction.
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$:

\[\begin{equation*} f^{-1}(t) = \sum_{n=1}^{N} \frac{t^n}{n!} \left. \frac{{\rm d}^{n-1}}{{\rm d} z^{n-1}}\left(\frac{z}{f(z)}\right)^n \right\vert_{z=0}. \end{equation*}\]
source
Base.absFunction.
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.

source
Base.LinAlg.normFunction.
norm(x::AbstractSeries, p::Real)

Returns the p-norm of an x::AbstractSeries, defined by

\[\begin{equation*} \left\Vert x \right\Vert_p = \left( \sum_k | x_k |^p \right)^{\frac{1}{p}}, \end{equation*}\]

which returns a non-negative number.

source
Base.isapproxFunction.
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.

source
Base.isfiniteFunction.
isfinite(x::AbstractSeries) -> Bool

Test whether the coefficients of the polynomial x are finite.

source
displayBigO(d::Bool) --> nothing

Set/unset displaying of the big 𝒪 notation in the output of Taylor1 and TaylorN polynomials. The initial value is true.

source

Internals

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.

source
generate_index_vectors(num_vars, degree)

Return a vector of index vectors with num_vars (number of variables) and degree.

source
TaylorSeries.in_baseFunction.
in_base(order, v)

Convert vector v of non-negative integers to base order+1.

source
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.

source
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.

source
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.

source
constant_term(a)

Return the constant value (zero order coefficient) for Taylor1 and TaylorN.

source
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

\[c_k = \sum_{j=0}^k a_j b_{k-j}.\]
source
TaylorSeries.mul!Method.
mul!(c, a, b) --> nothing

Return c = a*b with no allocation; all arguments are HomogeneousPolynomial.

source
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

\[c_k = \frac{1}{b_0} \big(a_k - \sum_{j=0}^{k-1} c_j b_{k-j}\big).\]

For Taylor1 polynomials, ordfact is the order of the factorized term of the denominator.

source
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

\[c_k = \frac{1}{k a_0} \sum_{j=0}^{k-1} \big(r(k-j) -j\big)a_{k-j} c_j.\]

For Taylor1 polynomials, k0 is the order of the first non-zero coefficient of a.

source
TaylorSeries.squareFunction.
square(a::AbstractSeries) --> typeof(a)

Return a^2; see TaylorSeries.sqr!.

source
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

\[\begin{eqnarray*} c_k & = & 2 \sum_{j=0}^{(k-1)/2} a_{k-j} a_j, \text{ if k is odd,} \\ c_k & = & 2 \sum_{j=0}^{(k-2)/2} a_{k-j} a_j + (a_{k/2})^2, \text{ if k is even. } \end{eqnarray*}\]
source
TaylorSeries.sqr!Method.
sqr!(c, a)

Return c = a*a with no allocation; all parameters are HomogeneousPolynomial.

source
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

\[\begin{eqnarray*} c_k &=& \frac{1}{2 c_0} \big( a_k - 2 \sum_{j=1}^{(k-1)/2} c_{k-j}c_j\big), \text{ if $k$ is odd,} \\ c_k &=& \frac{1}{2 c_0} \big( a_k - 2 \sum_{j=1}^{(k-2)/2} c_{k-j}c_j - (c_{k/2})^2\big), \text{ if $k$ is even.} \end{eqnarray*}\]

For Taylor1 polynomials, k0 is the order of the first non-zero coefficient, which must be even.

source
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

\[\begin{equation*} c_k = \frac{1}{k} \sum_{j=0}^{k-1} (k-j) a_{k-j} c_j. \end{equation*}\]
source
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

\[\begin{equation*} c_k = \frac{1}{a_0} \big(a_k - \frac{1}{k} \sum_{j=0}^{k-1} j a_{k-j} c_j \big). \end{equation*}\]
source
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

\[\begin{eqnarray*} s_k &=& \frac{1}{k}\sum_{j=0}^{k-1} (k-j) a_{k-j} c_j ,\\ c_k &=& -\frac{1}{k}\sum_{j=0}^{k-1} (k-j) a_{k-j} s_j. \end{eqnarray*}\]
source
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

\[\begin{equation*} c_k = a_k + \frac{1}{k} \sum_{j=0}^{k-1} (k-j) a_{k-j} p_j. \end{equation*}\]
source
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.

\[\begin{equation*} c_k = \frac{1}{ \sqrt{r_0} } \big( a_k - \frac{1}{k} \sum_{j=1}^{k-1} j r_{k-j} c_j \big). \end{equation*}\]
source
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.

\[\begin{equation*} c_k = - \frac{1}{ r_0 } \big( a_k - \frac{1}{k} \sum_{j=1}^{k-1} j r_{k-j} c_j \big). \end{equation*}\]
source
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.

\[\begin{equation*} c_k = \frac{1}{r_0}\big(a_k - \frac{1}{k} \sum_{j=1}^{k-1} j r_{k-j} c_j\big). \end{equation*}\]
source
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

\[\begin{eqnarray*} s_k &=& \frac{1}{k} \sum_{j=0}^{k-1} (k-j) a_{k-j} c_j, \\ c_k &=& \frac{1}{k} \sum_{j=0}^{k-1} (k-j) a_{k-j} s_j. \end{eqnarray*}\]
source
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.

\[\begin{equation*} c_k = a_k - \frac{1}{k} \sum_{j=0}^{k-1} (k-j) a_{k-j} p_j. \end{equation*}\]
source
Base.LinAlg.A_mul_B!Function.
A_mul_B!(Y, A, B)

Multiply A*B and save the result in Y.

source
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.

source
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

\[\begin{equation*} p_k = (k+1)a_{k+1}. \end{equation*}\]
source

Index