Library
Module
TaylorSeries.TaylorSeries
— ModuleTaylorSeries
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
— TypeTaylor1{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 :: Int
Maximum order (degree) of the polynomial.
Note that Taylor1
variables are callable. For more information, see evaluate
.
TaylorSeries.HomogeneousPolynomial
— TypeHomogeneousPolynomial{T<:Number} <: AbstractSeries{T}
DataType for homogeneous 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 homogeneous polynomial.
Note that HomogeneousPolynomial
variables are callable. For more information, see evaluate
.
TaylorSeries.TaylorN
— TypeTaylorN{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
— TypeAbstractSeries{T<:Number} <: Number
Parameterized abstract type for Taylor1
, HomogeneousPolynomial
and TaylorN
.
Functions and methods
TaylorSeries.Taylor1
— MethodTaylor1([T::Type=Float64], order::Int)
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
— MethodHomogeneousPolynomial([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
— MethodTaylorN([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
— Functionset_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{Int},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
— Functionget_variables(T::Type, [order::Int=get_order()])
Return a TaylorN{T}
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 explicitly established by the user without changing internal values for num_vars
or variable_names
. Omitting T
defaults to Float64
.
TaylorSeries.show_params_TaylorN
— Functionshow_params_TaylorN()
Display the current parameters for TaylorN
and HomogeneousPolynomial
types.
TaylorSeries.show_monomials
— Functionshow_monomials(ord::Int) --> nothing
List the indices and corresponding of a HomogeneousPolynomial
of degree ord
.
TaylorSeries.getcoeff
— Functiongetcoeff(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
, which is a tuple (or vector) with the indices of the specific monomial.
getcoeff(a, v)
Return the coefficient of a::TaylorN
, specified by v
, which is a tuple (or vector) with the indices of the specific monomial.
TaylorSeries.evaluate
— Functionevaluate(a, [dx])
Evaluate a Taylor1
polynomial using Horner's rule (hand coded). If dx
is omitted, its value is considered as zero. Note that 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::AbstractArray{Taylor1{T}}
, representing the dependent variables of an ODE, at time δt. Note that the syntax x(δt)
is equivalent to evaluate(x, δt)
, and x()
is equivalent to evaluate(x)
.
evaluate(a::Taylor1, x::Taylor1)
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 omitted, 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]; sorting::Bool=true)
Evaluate the TaylorN
polynomial a
at vals
. If vals
is omitted, it's evaluated at zero. The keyword parameter sorting
can be used to avoid sorting (in increasing order by abs2
) the terms that are added.
Note that the syntax a(vals)
is equivalent to evaluate(a, vals)
; and a()
is equivalent to evaluate(a)
; use a(b::Bool, x) corresponds to evaluate(a, x, sorting=b).
TaylorSeries.evaluate!
— Functionevaluate!(x, δt, x0)
Evaluates each element of x::AbstractArray{Taylor1{T}}
, 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
— Functiontaylor_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!
— Functionupdate!(a, x0)
Takes a <: Union{Taylo1,TaylorN}
and expands it around the coordinate x0
.
TaylorSeries.differentiate
— Functiondifferentiate(a)
Return the Taylor1
polynomial of the differential of a::Taylor1
. The order of the result is a.order-1
.
The function derivative
is an exact synonym of differentiate
.
differentiate(a, n)
Compute recursively the Taylor1
polynomial of the n-th derivative of a::Taylor1
. The order of the result is a.order-n
.
differentiate(n, a)
Return the value of the n
-th differentiate of the polynomial a
.
differentiate(a, r)
Partial differentiation of a::HomogeneousPolynomial
series with respect to the r
-th variable.
differentiate(a, r)
Partial differentiation of a::TaylorN
series with respect to the r
-th variable. The r
-th variable may be also specified through its symbol.
differentiate(a::TaylorN{T}, ntup::NTuple{N,Int})
Return a TaylorN
with the partial derivative of a
defined by ntup::NTuple{N,Int}
, where the first entry is the number of derivatives with respect to the first variable, the second is the number of derivatives with respect to the second, and so on.
differentiate(ntup::NTuple{N,Int}, a::TaylorN{T})
Returns the value of the coefficient of a
specified by ntup::NTuple{N,Int}
, multiplied by the corresponding factorials.
TaylorSeries.derivative
— Functionderivative
An exact synonym of differentiate
.
TaylorSeries.integrate
— Functionintegrate(a, [x])
Return the integral of a::Taylor1
. The constant of integration (0-th order coefficient) is set to x
, which is zero if omitted. Note that the order of the result is a.order+1
.
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 omitted, it is taken as zero.
TaylorSeries.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 omitted, 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 omitted, 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 omitted, 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 omitted, it is evaluated at zero.
TaylorSeries.constant_term
— Functionconstant_term(a)
Return the constant value (zero order coefficient) for Taylor1
and TaylorN
. The fallback behavior is to return a
itself if a::Number
, or a[1]
when a::Vector
.
TaylorSeries.linear_polynomial
— Functionlinear_polynomial(a)
Returns the linear part of a
as a polynomial (Taylor1
or TaylorN
), without the constant term. The fallback behavior is to return a
itself.
TaylorSeries.nonlinear_polynomial
— Functionnonlinear_polynomial(a)
Returns the nonlinear part of a
. The fallback behavior is to return zero(a)
.
TaylorSeries.inverse
— Functioninverse(f)
Return the Taylor expansion of $f^{-1}(t)$, of order N = f.order
, for f::Taylor1
polynomial, assuming the first coefficient of f
is zero. Otherwise, a DomainError
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*}\]
TaylorSeries.inverse_map
— Functioninverse_map(f)
Return the Taylor expansion of $f^{-1}(t)$, of order N = f.order
, for Taylor1
or TaylorN
polynomials, assuming the first coefficient of f
is zero. Otherwise, a DomainError
is thrown.
This method is based in the algorithm by M. Berz, Modern map methods in Particle Beam Physics, Academic Press (1999), Sect 2.3.1. See inverse
(for f::Taylor1
).
Base.abs
— Functionabs(a)
For a Real
type returns a
if constant_term(a) > 0
and -a
if constant_term(a) < 0
for a <:Union{Taylor1,TaylorN}
. For a Complex
type, such as Taylor1{ComplexF64}
, returns sqrt(real(a)^2 + imag(a)^2)
.
Notice that typeof(abs(a)) <: AbstractSeries
and that for a Complex
argument a Real
type is returned (e.g. typeof(abs(a::Taylor1{ComplexF64})) == Taylor1{Float64}
).
LinearAlgebra.norm
— Functionnorm(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.
Base.isapprox
— Functionisapprox(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.isless
— Functionisless(a::Taylor1{<:Real}, b::Real)
isless(a::TaylorN{<:Real}, b::Real)
Compute isless
by comparing the constant_term(a)
and b
. If they are equal, returns a[nz] < 0
, with nz
the first non-zero coefficient after the constant term. This defines a total order.
For many variables, the ordering includes a lexicographical convention in order to be total. We have opted for the simplest one, where the larger variable appears before when the TaylorN
variables are defined (e.g., through set_variables
).
Refs:
- M. Berz, AIP Conference Proceedings 177, 275 (1988); https://doi.org/10.1063/1.37800
- M. Berz, "Automatic Differentiation as Nonarchimedean Analysis", Computer Arithmetic and Enclosure Methods, (1992), Elsevier, 439-450.
isless(a::Taylor1{<:Real}, b::Taylor1{<:Real})
isless(a::TaylorN{<:Real}, b::Taylor1{<:Real})
Returns isless(a - b, zero(b))
.
Base.isfinite
— Functionisfinite(x::AbstractSeries) -> Bool
Test whether the coefficients of the polynomial x
are finite.
TaylorSeries.displayBigO
— FunctiondisplayBigO(d::Bool) --> nothing
Set/unset displaying of the big 𝒪 notation in the output of Taylor1
and TaylorN
polynomials. The initial value is true
.
TaylorSeries.use_show_default
— Functionuse_Base_show(d::Bool) --> nothing
Use Base.show_default
method (default show
method in Base), or a custom display. The initial value is false
, so customized display is used.
TaylorSeries.set_taylor1_varname
— Functionset_taylor1_varname(var::String)
Change the displayed variable for Taylor1
objects.
Internals
TaylorSeries.ParamsTaylor1
— TypeParamsTaylor1
DataType holding the current variable name for Taylor1
.
Field:
var_name :: String
Names of the variables
These parameters can be changed using set_taylor1_varname
TaylorSeries.ParamsTaylorN
— TypeParamsTaylorN
DataType holding the current parameters for TaylorN
and HomogeneousPolynomial
.
Fields:
order :: Int
Order (degree) of the polynomialsnum_vars :: Int
Number of variablesvariable_names :: Vector{String}
Names of the variablesvariable_symbols :: Vector{Symbol}
Symbols of the variables
These parameters can be changed using set_variables
TaylorSeries._InternalMutFuncs
— Type_InternalMutFuncs
Contains parameters and expressions that allow a simple programmatic construction for calling the internal mutating functions.
TaylorSeries.generate_tables
— Functiongenerate_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{Int,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{Int,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{Int,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{Int,Int},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
— Functiongenerate_index_vectors(num_vars, degree)
Return a vector of index vectors with num_vars
(number of variables) and degree.
TaylorSeries.in_base
— Functionin_base(order, v)
Convert vector v
of non-negative integers to base oorder
, where oorder
is the next odd integer of order
.
TaylorSeries.make_inverse_dict
— Functionmake_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!
— Functionresize_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!
— Functionresize_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.numtype
— Functionnumtype(a::AbstractSeries)
Returns the type of the elements of the coefficients of a
.
LinearAlgebra.mul!
— Functionmul!(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
. Note that for TaylorN
the result of a * b
is accumulated in c[k]
.
The coefficients are given by
\[c_k = \sum_{j=0}^k a_j b_{k-j}.\]
LinearAlgebra.mul!
— Methodmul!(c, a, b) --> nothing
Accumulates in c
the result of a*b
with minimum allocation. Arguments c, a and b are HomogeneousPolynomial
.
TaylorSeries.mul_scalar!
— Methodmul_scalar!(c, scalar, a, b) --> nothing
Accumulates in c
the result of scalar*a*b
with minimum allocation. Arguments c, a and b are HomogeneousPolynomial
; scalar
is a NumberNotSeries.
LinearAlgebra.mul!
— Methodmul!(Y, A, B)
Multiply A*B and save the result in Y.
TaylorSeries.div!
— Functiondiv!(c, a, b, k::Int)
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, a similar formula is implemented which exploits k_0
, the order of the first non-zero coefficient of a
.
TaylorSeries.pow!
— Functionpow!(c, a, r::Real, k::Int)
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, a similar formula is implemented which exploits k_0
, the order of the first non-zero coefficient of a
.
TaylorSeries.square
— Functionsquare(a::AbstractSeries) --> typeof(a)
Return a^2
; see TaylorSeries.sqr!
.
TaylorSeries.sqr!
— Functionsqr!(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{aligned} 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{aligned}\]
TaylorSeries.accsqr!
— Functionaccsqr!(c, a)
Returns c += a*a
with no allocation; all parameters are HomogeneousPolynomial
.
TaylorSeries.sqrt!
— Functionsqrt!(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{aligned} 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{aligned}\]
For Taylor1
polynomials, k0
is the order of the first non-zero coefficient, which must be even.
TaylorSeries.exp!
— Functionexp!(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*}\]
TaylorSeries.log!
— Functionlog!(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*}\]
TaylorSeries.sincos!
— Functionsincos!(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{aligned} 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{aligned}\]
TaylorSeries.tan!
— Functiontan!(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*}\]
TaylorSeries.asin!
— Functionasin!(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*}\]
TaylorSeries.acos!
— Functionacos!(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*}\]
TaylorSeries.atan!
— Functionatan!(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*}\]
TaylorSeries.sinhcosh!
— Functionsinhcosh!(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{aligned} 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{aligned}\]
TaylorSeries.tanh!
— Functiontanh!(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*}\]
TaylorSeries.asinh!
— Functionasinh!(c, a, r, k)
Update the k-th
expansion coefficients c[k+1]
of c = asinh(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*}\]
TaylorSeries.acosh!
— Functionacosh!(c, a, r, k)
Update the k-th
expansion coefficients c[k+1]
of c = acosh(a)
, for c
and a
either Taylor1
or TaylorN
; r = sqrt(c^2-1)
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*}\]
TaylorSeries.atanh!
— Functionatanh!(c, a, r, k)
Update the k-th
expansion coefficients c[k+1]
of c = atanh(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*}\]
TaylorSeries.differentiate!
— Functiondifferentiate!(res, a) --> nothing
In-place version of differentiate
. Compute the Taylor1
polynomial of the differential of a::Taylor1
and return it as res
(order of res
remains unchanged).
differentiate!(p, a, k) --> nothing
Update in-place the k-th
expansion coefficient p[k]
of p = differentiate(a)
for both p
and a
Taylor1
.
The coefficients are given by
\[p_k = (k+1) a_{k+1}.\]
TaylorSeries._internalmutfunc_call
— Function_internalmutfunc_call( fn :: _InternalMutFuncs )
Creates the appropriate call to the internal mutating function defined by the _InternalMutFuncs
object. This is used to construct _dict_unary_calls
and _dict_binary_calls
. The call contains the prefix TaylorSeries.
.
TaylorSeries._dict_unary_ops
— Constant_dict_binary_ops
Dict{Symbol, Array{Any,1}}
with the information to construct the _InternalMutFuncs
related to unary operations.
The keys correspond to the function symbols.
The arguments of the array are the function name (e.g. add!
), a tuple with the function arguments, and an Expr
with the calling pattern. The convention for the arguments of the functions and the calling pattern is to use :_res
for the (mutated) result, :_arg1
, for the required argument, possibly :_aux
when there is an auxiliary expression needed, and :_k
for the computed order of :_res
. When an auxiliary expression is required, an Expr
defining its calling pattern is added as the last entry of the vector.
TaylorSeries._dict_binary_calls
— Constant_dict_binary_calls::Dict{Symbol, NTuple{2,Expr}}
Dictionary with the expressions that define the internal binary functions and the auxiliary functions, whenever they exist. The keys correspond to those functions, passed as symbols, with the defined internal mutating functions.
Evaluating the entries generates symbols that represent the actual calls to the internal mutating functions.
TaylorSeries._dict_unary_calls
— Constant_dict_unary_calls::Dict{Symbol, NTuple{2,Expr}}
Dictionary with the expressions that define the internal unary functions and the auxiliary functions, whenever they exist. The keys correspond to those functions, passed as symbols, with the defined internal mutating functions.
Evaluating the entries generates expressions that represent the actual calls to the internal mutating functions.
TaylorSeries._dict_binary_ops
— Constant_dict_binary_ops
Dict{Symbol, Array{Any,1}}
with the information to construct the _InternalMutFuncs
related to binary operations.
The keys correspond to the function symbols.
The arguments of the array are the function name (e.g. add!
), a tuple with the function arguments, and an Expr
with the calling pattern. The convention for the arguments of the functions and the calling pattern is to use :_res
for the (mutated) result, :_arg1
and _arg2
for the required arguments, and :_k
for the computed order of :_res
.
TaylorSeries._populate_dicts!
— Function_populate_dicts!()
Function that populates the internal dictionaries _dict_unary_calls
and _dict_binary_calls
TaylorSeries.@isonethread
— Macro@isonethread (expr)
Internal macro used to check the number of threads in use, to prevent a data race that modifies coeff_table
when using differentiate
or integrate
; see https://github.com/JuliaDiff/TaylorSeries.jl/issues/318.
This macro is inspired by the macro @threaded
; see https://github.com/trixi-framework/Trixi.jl/blob/main/src/auxiliary/auxiliary.jl; and https://github.com/trixi-framework/Trixi.jl/pull/426/files.
Index
TaylorSeries.AbstractSeries
TaylorSeries.HomogeneousPolynomial
TaylorSeries.HomogeneousPolynomial
TaylorSeries.ParamsTaylor1
TaylorSeries.ParamsTaylorN
TaylorSeries.Taylor1
TaylorSeries.Taylor1
TaylorSeries.TaylorN
TaylorSeries.TaylorN
TaylorSeries._InternalMutFuncs
Base.abs
Base.isapprox
Base.isfinite
Base.isless
LinearAlgebra.mul!
LinearAlgebra.mul!
LinearAlgebra.mul!
LinearAlgebra.norm
TaylorSeries._internalmutfunc_call
TaylorSeries._populate_dicts!
TaylorSeries.accsqr!
TaylorSeries.acos!
TaylorSeries.acosh!
TaylorSeries.asin!
TaylorSeries.asinh!
TaylorSeries.atan!
TaylorSeries.atanh!
TaylorSeries.constant_term
TaylorSeries.derivative
TaylorSeries.differentiate
TaylorSeries.differentiate!
TaylorSeries.displayBigO
TaylorSeries.div!
TaylorSeries.evaluate
TaylorSeries.evaluate!
TaylorSeries.exp!
TaylorSeries.generate_index_vectors
TaylorSeries.generate_tables
TaylorSeries.get_variables
TaylorSeries.getcoeff
TaylorSeries.gradient
TaylorSeries.hessian
TaylorSeries.hessian!
TaylorSeries.in_base
TaylorSeries.integrate
TaylorSeries.inverse
TaylorSeries.inverse_map
TaylorSeries.jacobian
TaylorSeries.jacobian!
TaylorSeries.linear_polynomial
TaylorSeries.log!
TaylorSeries.make_inverse_dict
TaylorSeries.mul_scalar!
TaylorSeries.nonlinear_polynomial
TaylorSeries.numtype
TaylorSeries.pow!
TaylorSeries.resize_coeffs1!
TaylorSeries.resize_coeffsHP!
TaylorSeries.set_taylor1_varname
TaylorSeries.set_variables
TaylorSeries.show_monomials
TaylorSeries.show_params_TaylorN
TaylorSeries.sincos!
TaylorSeries.sinhcosh!
TaylorSeries.sqr!
TaylorSeries.sqrt!
TaylorSeries.square
TaylorSeries.tan!
TaylorSeries.tanh!
TaylorSeries.taylor_expand
TaylorSeries.update!
TaylorSeries.use_show_default