Library¶
Types¶
#
TaylorSeries.Taylor1
— Type.
Taylor1{T<:Number} <: Number
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.
#
TaylorSeries.HomogeneousPolynomial
— Type.
HomogeneousPolynomial{T<:Number} <: Number
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.
#
TaylorSeries.TaylorN
— Type.
TaylorN{T<:Number} <: Number
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.
#
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, [order=1])
Taylor1([order=1])
Shortcut to define the independent variable of a Taylor1{T}
polynomial of given order
. If the type T
is ommitted, Float64
is assumend.
#
TaylorSeries.TaylorN
— Method.
TaylorN(T, nv; [order=get_order()])
TaylorN(nv; [order=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()
. If T::Type
is ommitted, Float64
is assumend.
#
TaylorSeries.set_variables
— Function.
set_variables(R, names; [order=6])
set_variables(names; [order=6])
set_variables(R, names; [order=6, numvars=-1])
set_variables(names; [order=6, numvars=-1])
Return a TaylorN{R}
vector with numvars
polynomials, each representing an independent variable, using names
as the String
for the output.
If numvars
is not specified, it is inferred from the length of names
. If length(names)==1
and numvars>1
, it uses this name with subscripts for the different variables. When changing the order
or numvars
, the hash_tables are reset.
#
TaylorSeries.show_params_TaylorN
— Function.
show_params_TaylorN()
Display the current parameters for TaylorN
and HomogeneousPolynomial
types.
#
TaylorSeries.get_coeff
— Function.
get_coeff(a, n)
Return the coefficient of order n::Int
of a a::Taylor1
polynomial.
get_coeff(a, v)
Return the coefficient of a::HomogeneousPolynomial
, specified by v::Array{Int,1}
which has the indices of the specific monomial.
get_coeff(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.
evaluate(x, δt)
Evaluates each element of x::Array{Taylor1{T},1}
, representing the dependent variables of an ODE, at time δt.
evaluate(a, x)
Substitute x::Taylor1
as independent variable in a a::Taylor1
polynomial.
evaluate(a, vals)
Evaluate a HomogeneousPolynomial
polynomial using Horner's rule (hand coded) at vals
.
evaluate(a, [vals])
Evaluate a TaylorN
polynomial using Horner's rule (hand coded) at vals
. If vals
is ommitted, it is evaluated at zero.
#
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.derivative
— Function.
derivative(a)
Return the Taylor1
polynomial of the differential of a::Taylor1
. The last coefficient is set to zero.
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.
#
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.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.
#
Base.:*
— Function.
*(a, b)
Return the Taylor expansion of $a \cdot b$, of order max(a.order,b.order)
, for a::Taylor1
, b::Taylor1
polynomials.
For details on making the Taylor expansion, see TaylorSeries.mulHomogCoef
.
#
Base.:/
— Function.
/(a, b)
Return the Taylor expansion of $a/b$, of order max(a.order,b.order)
, for a::Taylor1
, b::Taylor1
polynomials.
For details on making the Taylor expansion, see TaylorSeries.divHomogCoef
.
#
Base.:^
— Function.
^(a, x)
Return the Taylor expansion of $a^x$ for a::Taylor1
polynomial and x::Number
. If x
is non integer and the 0-th order coefficient is zero, an ArgumentError
is thrown.
#
Base.sqrt
— Function.
sqrt(a)
Return the Taylor expansion of $\sqrt(a)$, of order a.order
, for a::Taylor1
polynomial. If the first non-vanishing coefficient of a
is an odd power, and ArgumentError
is thrown.
For details on making the Taylor expansion, see TaylorSeries.sqrtHomogCoef
.
#
Base.exp
— Function.
exp(a)
Return the Taylor expansion of $e^a$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.expHomogCoef
.
#
Base.log
— Function.
log(a)
Return the Taylor expansion of $\log(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.logHomogCoef
.
#
Base.sin
— Function.
sin(a)
Return the Taylor expansion of $\sin(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.sincosHomogCoef
.
#
Base.cos
— Function.
cos(a)
Return the Taylor expansion of $\cos(a)$, of order a.order
, for a::Taylor1
polynomial
For details on making the Taylor expansion, see TaylorSeries.sincosHomogCoef
.
#
Base.tan
— Function.
tan(a)
Return the Taylor expansion of $\tan(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.tanHomogCoef
.
#
Base.asin
— Function.
asin(a)
Return the Taylor expansion of $\arcsin(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.asinHomogCoef
.
#
Base.acos
— Function.
acos(a)
Return the Taylor expansion of $\arccos(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.acosHomogCoef
.
#
Base.atan
— Function.
atan(a)
Return the Taylor expansion of $\arctan(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.atanHomogCoef
.
#
Base.sinh
— Function.
sinh(a)
Return the Taylor expansion of $\sinh(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.sinhcoshHomogCoef
.
#
Base.cosh
— Function.
cosh(a)
Return the Taylor expansion of $\cosh(a)$, of order a.order
, for a::Taylor1
polynomial
For details on making the Taylor expansion, see TaylorSeries.sinhcoshHomogCoef
.
#
Base.tanh
— Function.
tanh(a)
Return the Taylor expansion of $\tanh(a)$, of order a.order
, for a::Taylor1
polynomial.
For details on making the Taylor expansion, see TaylorSeries.tanhHomogCoef
.
#
Base.abs
— Function.
abs(a)
Return a
or -a
depending on the 0-th order coefficient of the Taylor1
polynomial a
. If a.coeffs[1]
is zero, an ArgumentError
is thrown.
abs(a::TaylorN)
Absolute value of a TaylorN
polynomial, using the 0-th order coefficient.
Return a
or -a
, depending on the 0-th order coefficient of a
. If it is zero, it throws an ArgumentError
.
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.order_posTb
— Function.
order_posTb(order, nv, ord)
Return a vector with the positions, in a HomogeneousPolynomial
of order order
, where the variable nv
has order ord
.
#
TaylorSeries.mul!
— Function.
mul!(c, a, b)
Return c = a*b
with no allocation; all parameters are HomogeneousPolynomial
.
mul!(c, a)
Return c = a*a
with no allocation; all parameters are HomogeneousPolynomial
.
#
TaylorSeries.mulHomogCoef
— Function.
mulHomogCoef(kcoef, ac, bc)
Compute the k
-th expansion coefficient of $c = a\cdot b$ given by
with $a$ and $b$ Taylor1
polynomials.
Inputs are the kcoef
-th coefficient, and the vectors of the expansion coefficients ac
and bc
, corresponding respectively to a
and b
.
#
TaylorSeries.divHomogCoef
— Function.
divHomogCoef(kcoef, ac, bc, coeffs, ordfact)
Compute the k-th
expansion coefficient of $c = a / b$ given by
with $a$ and $b$ Taylor1
polynomials.
Inputs are the kcoef
-th coefficient, the vectors of the expansion coefficients ac
and bc
, corresponding respectively to a
and b
, the already calculated expansion coefficients coeffs
of c
, and ordfact
which is the order of the factorized term of the denominator, whenever b_0
is zero.
#
TaylorSeries.powHomogCoef
— Function.
powHomogCoef(kcoef, ac, x, coeffs, knull)
Compute the k-th
expansion coefficient of $c = a^x$, given by
with $a$ a Taylor1
polynomial, and x
a number.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, the exponent x
, the already calculated expansion coefficients coeffs
of c
, and knull
which is the order of the first non-zero coefficient of a
.
#
TaylorSeries.squareHomogCoef
— Function.
squareHomogCoef(kcoef, ac)
Compute the k-th
expansion coefficient of $c = a^2$, given by
with $a$ a Taylor1
polynomial.
Inputs are the kcoef
-th coefficient and the vector of the expansion coefficients ac
of a
.
#
TaylorSeries.sqrtHomogCoef
— Function.
sqrtHomogCoef(kcoef, ac, coeffs, knull)
Compute the k-th
expansion coefficient of $c = \sqrt(a)$, given by
with $a$ a Taylor1
polynomial.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, the already calculated expansion coefficients coeffs
of c
, and knull
, which is half of the order of the first non-zero coefficient of a
.
#
TaylorSeries.expHomogCoef
— Function.
expHomogCoef(kcoef, ac, coeffs)
Compute the k-th
expansion coefficient of $c = \exp(a)$ given by
with $a$ a Taylor1
polynomial.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of $a$, and the already calculated expansion coefficients coeffs
of c
.
#
TaylorSeries.logHomogCoef
— Function.
logHomogCoef(kcoef, ac, coeffs)
Compute the k-th
expansion coefficient of $c = \log(a)$, given by
with $a$ a Taylor1
polynomial.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, and the already calculated expansion coefficients coeffs
of c
.
#
TaylorSeries.sincosHomogCoef
— Function.
sincosHomogCoef(kcoef, ac, scoeffs, ccoeffs)
Compute the k-th
expansion coefficient of $s = \sin(a)$ and $c=\cos(a)$ simultaneously given by
with $a$ a Taylor1
polynomial.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, and the already calculated expansion coefficients scoeffs
and ccoeffs
of sin(a)
and cos(a)
, respectvely.
#
TaylorSeries.tanHomogCoef
— Function.
tanHomogCoef(kcoef, ac, coeffst2)
Compute the k-th
expansion coefficient of $c = \tan(a)$ given by
with $a$ a Taylor1
polynomial and $p = c^2$.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, and the already calculated expansion coefficients coeffst2
of c^2
.
#
TaylorSeries.asinHomogCoef
— Function.
asinHomogCoef(kcoef, ac, rc, coeffs)
Compute the k-th
expansion coefficient of $s = \arcsin(a)$ given by
with $a$ a Taylor1
polynomial and $r = \sqrt{1 - a^2}$.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, the already calculated expansion coefficients rc
of $r$ (see above), and the already calculated expansion coefficients coeffs
of asin(a)
.
#
TaylorSeries.acosHomogCoef
— Function.
acosHomogCoef(kcoef, ac, rc, coeffs)
Compute the k-th
expansion coefficient of $c = \arccos(a)$ given by
with $a$ a Taylor1
polynomial and $r = \sqrt{1 - a^2}$.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, the already calculated expansion coefficients rc
of $r$ (see above), and the already calculated expansion coefficients coeffs
of acos(a)
.
#
TaylorSeries.atanHomogCoef
— Function.
atanHomogCoef(kcoef, ac, rc, coeffs)
Compute the k-th
expansion coefficient of $c = \arctan(a)$ given by
with $a$ a Taylor1
polynomial and $r = 1 + a^2$.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, the already calculated expansion coefficients rc
of $r$ (see above), and the already calculated expansion coefficients coeffs
of asin(a)
.
#
TaylorSeries.sinhcoshHomogCoef
— Function.
sinhcoshHomogCoef(kcoef, ac, scoeffs, ccoeffs)
Compute the k-th
expansion coefficient of $sh = \sinh(a)$ and $ch=\cosh(a)$ simultaneously given by
with $a$ a Taylor1
polynomial.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, and the already calculated expansion coefficients shcoeffs
and chcoeffs
of sinh(a)
and cosh(a)
, respectvely.
#
TaylorSeries.tanhHomogCoef
— Function.
tanhHomogCoef(kcoef, ac, coeffst2)
Compute the k-th
expansion coefficient of $c = \tanh(a)$ given by
with $a$ a Taylor1
polynomial and $p = th^2$.
Inputs are the kcoef
-th coefficient, the vector of the expansion coefficients ac
of a
, and the already calculated expansion coefficients coeffst2
of th^2
.
#
Base.LinAlg.A_mul_B!
— Function.
A_mul_B!(Y, A, B)
Multiply A*B and save the result in Y.