GTPSA.CTPSAType
`CTPSA`

This is a 1-to-1 struct for the C definition ctpsa (complex TPSA) in GTPSA.

Fields

  • d::Ptr{Cvoid} – Pointer to Desc for this CTPSA
  • uid::Cint – Special user field for external use (and padding)
  • mo::Cuchar – max ord (allocated)
  • lo::Cuchar – lowest used ord
  • hi::Cuchar – highest used ord
  • nz::Culonglong – zero/nonzero homogenous polynomials. Note: Int64 if 64 bit compiled C code, else 32 bit
  • nam::NTuple{NAMSZ,Cuchar} – tpsa name max string length 16 NAMSZ
  • coef::Ptr{ComplexF64} – warning: must be identical to ctpsa up to coef excluded
GTPSA.ComplexTPSType
ComplexTPS(cta::Union{Number,Nothing}=nothing; use::Union{Descriptor,TPS,ComplexTPS,Nothing}=nothing)::ComplexTPS

Constructor to create a new ComplexTPS equal to the number cta. If cta is a ComplexTPS (or TPS), this is equivalent to a copy constructor, with the result by default having the same Descriptor as cta. If cta is not a TPS orComplexTPS, then the Descriptor used will by default be GTPSA.desc_current. The Descriptor for the constructed ComplexTPS can be set using use. If a TPS or ComplexTPS is passed to use, the Descriptor for that TPS will be used.

The constructor can also be used to create a copy of a ComplexTPS under one Descriptor to instead have a different Descriptor. In this case, invalid monomials under the new Descriptor are removed.

Input

  • cta – Any Number
  • use – (Optional) specify which Descriptor to use, default is nothing which uses the Descriptor for cta if cta <: Union{TPS,ComplexTPS}, else uses GTPSA.desc_current

Output

  • ret – New ComplexTPS equal to cta with removal of invalid monomials if cta is a TPS/ComplexTPS and a new Descriptor is specified
GTPSA.DescType
`Desc`

This is a 1-to-1 struct for the C definition desc (descriptor) in GTPSA. Descriptors include all information about the TPSA, including the number of variables/parameters and their orders, lookup tables for the monomials, monomial indexing function, and pre-allocated permanent temporaries for fast evaluation.

Fields

  • id::Cint – Index in list of registered descriptors
  • nn::Cint – Number of variables + number of parameters, nn = nv+np <= 100000
  • nv::Cint – Number of variables
  • np::Cint – Number of parameters
  • mo::Cuchar – Max order of both variables AND parameters
  • po::Cuchar – Max order of parameterss
  • to::Cuchar – Global order of truncation. Note: ord_t in gtpsa is typedef for unsigned char (Cuchar)
  • no::Ptr{Cuchar} – Array of orders of each variable (first nv entries) and parameters (last np entries), length nn. Note: In C this is const
  • uno::Cint – User provided array of orders of each variable/parameter (with mad_desc_newvpo)
  • nth::Cint – Max number of threads or 1
  • nc::Cuint – Number of coefficients (max length of TPSA)
  • monos::Ptr{Cuchar} – 'Matrix' storing the monomials (sorted by variable)
  • ords::Ptr{Cuchar} – Order of each monomial of To
  • prms::Ptr{Cuchar} – Order of parameters in each monomial of To (zero = no parameters)
  • To::Ptr{Ptr{Cuchar}} – Table by orders - pointers to monomials, sorted by order
  • Tv::Ptr{Ptr{Cuchar}} – Table by vars - pointers to monomials, sorted by variable
  • ocs::Ptr{Ptr{Cuchar}}ocs[t,i] -> o in mul, compute o on thread t 3 <= o <= mo aterminated with 0
  • ord2idx::Ptr{Cint} – Order to polynomial start index in To (i.e. in TPSA coef)
  • tv2to::Ptr{Cint} – Lookup tv->to
  • to2tv::Ptr{Cint} – Lookup to->tv
  • H::Ptr{Cint} – Indexing matrix in Tv
  • L::Ptr{Ptr{Cint}} – Multiplication indexes L[oa,ob]->L_ord L_ord[ia,ib]->ic
  • L_idx::Ptr{Ptr{Ptr{Cint}}}L_idx[oa,ob]->[start] [split] [end] idxs in L
  • size::Culonglong – Bytes used by desc. Unsigned Long Int: In 32 bit system is Int32 but 64 bit is Int64. Using Culonglong assuming 64 bit
  • t::Ptr{Ptr{Cvoid}} – Temporary array contains 8 pointers to RTPSAs already initialized
  • ct::Ptr{Ptr{Cvoid}} – Temporary array contains 8 pointers to CTPSAs already initialized
  • ti::Ptr{Cint} – idx of tmp used by each thread (length = # threads)
  • cti::Ptr{Cint} – idx of tmp used by each thread (length = # threads)
GTPSA.DescriptorMethod
Descriptor(nv::Integer, vo::Integer, np::Integer, po::Integer)::Descriptor

Creates a TPSA Descriptor with nv variables each with truncation order vo, and np parameters each with truncation order po

Input

  • nv – Number of variables in the TPSA
  • vo – Truncation order of the variables in the TPSA
  • np – Number of parameters in the TPSA
  • po – Truncation order of the parameters
GTPSA.DescriptorMethod
Descriptor(nv::Integer, vo::Integer)::Descriptor

Creates a TPSA Descriptor with nv variables of maximum order vo for each.

Input

  • nv – Number of variables in the TPSA
  • vo – Maximum order of the variables in the TPSA
GTPSA.DescriptorMethod
Descriptor(vos::Vector{<:Integer}, pos::Vector{<:Integer})::Descriptor

Creates a TPSA Descriptor with length(vos) variables with individual truncation orders specified in vos, and length(pos) parameters with individual truncation orders specified in pos.

Input

  • vosVector of the individual truncation orders of each variable
  • posVector of the individual truncation orders of each parameter
GTPSA.DescriptorMethod
Descriptor(vos::Vector{<:Integer})::Descriptor

Creates a TPSA Descriptor with length(mos) variables with individual truncation orders specified in the Vector vos.

Input

  • vosVector of the individual truncation orders of each variable
GTPSA.RTPSAType
`RTPSA`

This is a 1-to-1 struct for the C definition tpsa (real TPSA) in GTPSA.

Fields

  • d::Ptr{Cvoid} – Ptr to tpsa descriptor
  • uid::Cint – Special user field for external use (and padding)
  • mo::Cuchar – max ord (allocated)
  • lo::Cuchar – lowest used ord
  • hi::Cuchar – highest used ord
  • nz::Culonglong – zero/nonzero homogenous polynomials. Note: Int64 if 64 bit compiled C code, else 32 bit
  • nam::NTuple{NAMSZ,Cuchar} – tpsa name max string length 16 NAMSZ
  • coef::Ptr{Cdouble} – warning: must be identical to ctpsa up to coef excluded
GTPSA.TPSType
TPS(ta::Union{Real,Nothing}=nothing; use::Union{Descriptor,TPS,ComplexTPS,Nothing}=nothing)::TPS

Constructor to create a new TPS equal to the real value ta. If ta is a TPS, this is equivalent to a copy constructor, with the result by default having the same Descriptor as ta. If ta is not a TPS, then the Descriptor used will by default be GTPSA.desc_current. The Descriptor for the constructed TPS can be set using use. If a TPS or ComplexTPS is passed to use, the Descriptor for that TPS will be used.

The constructor can also be used to create a copy of a TPS under one Descriptor to instead have a different Descriptor. In this case, invalid monomials under the new Descriptor are removed.

Input

  • ta – Any Real
  • use – (Optional) specify which Descriptor to use, default is nothing which uses the Descriptor for ta if ta isa TPS, else uses GTPSA.desc_current

Output

  • ret – New TPS equal to ta with removal of invalid monomials if ta is a TPS and a new Descriptor is specified
Base.invMethod
inv(ma::Vector{<:Union{TPS,ComplexTPS}})

Inverts the map ma such that ma ∘ inv(ma) = 1 in the variables.

Example


julia> d = Descriptor(2,10); x = vars()[1]; p = vars()[2];

julia> time = 0.01; k = 2; m = 0.01;

julia> h = p^2/(2m) + 1/2*k*x^2;

julia> hf = getvectorfield(h);

julia> map = exppb(-time*hf, [x, p]);

julia> map ∘ inv(map)
2-element Vector{TPS}:
  Out  Coefficient                Order   Exponent
-------------------------------------------------
   1:   1.0000000000000000e+00      1      1   0
-------------------------------------------------
   2:   1.0000000000000002e+00      1      0   1
GTPSA.complexparamsFunction
complexparams(d::Descriptor=GTPSA.desc_current)::Vector{ComplexTPS}

Returns ComplexTPSs corresponding to the parameters for the Descriptor specified by use. Default value is GTPSA.desc_current.

Input

  • use – (Optional) Specify which TPSA Descriptor to use, default is GTPSA.desc_current

Output

  • kVector containing unit ComplexTPSs corresponding to each parameter
GTPSA.complexvarsFunction
complexvars(use::Union{Descriptor,TPS,ComplexTPS}=GTPSA.desc_current)::Vector{ComplexTPS}

Returns ComplexTPSs corresponding to the variables for the Descriptor specified by use. Default value is GTPSA.desc_current.

Input

  • use – (Optional) Specify which TPSA Descriptor to use, default is GTPSA.desc_current

Output

  • xVector containing unit ComplexTPSs corresponding to each variable
GTPSA.cutordMethod
cutord(t1::Union{TPS, ComplexTPS}, order::Integer)

Cuts out the monomials in t1 at the given order and above. Or, if order is negative, will cut monomials with orders at and below abs(order).

Examples

julia> d = Descriptor(1,10);

julia> x = vars(d);

julia> cutord(sin(x[1]), 5)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1
  -1.6666666666666666e-01    3        3


julia> cutord(sin(x[1]), -5)
TPS:
  Coefficient              Order     Exponent
  -1.9841269841269841e-04    7        7
   2.7557319223985893e-06    9        9
GTPSA.derivFunction
deriv(t1::Union{TPS,ComplexTPS}, v::Union{Integer, Vector{<:Union{<:Pair{<:Integer,<:Integer},<:Integer}}, Nothing}=nothing; param::Union{<:Integer,Nothing}=nothing, params::Union{Vector{<:Pair{<:Integer,<:Integer}}, Nothing}=nothing)
∂(t1::Union{TPS,ComplexTPS}, v::Union{Integer, Vector{<:Union{<:Pair{<:Integer,<:Integer},<:Integer}}, Nothing}=nothing; param::Union{<:Integer,Nothing}=nothing, params::Union{Vector{<:Pair{<:Integer,<:Integer}}, Nothing}=nothing)

Differentiates t1 wrt the variable/parameter specified by the variable/parameter index, or alternatively any monomial specified by indexing-by-order OR indexing-by-sparse monomial.

Input

  • v – An integer (for variable index), an array of orders for each variable (for indexing-by-order), or an array of pairs (sparse monomial)
  • param – (Keyword argument, optional) An integer for the parameter index
  • params – (Keyword argument, optional) An array of pairs for sparse-monomial indexing

Examples: Variable/Parameter Index:

julia> d = Descriptor(1,5,1,5);

julia> x1 = vars(d)[1]; k1 = params(d)[1];

julia> deriv(x1*k1, 1)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    1


julia> deriv(x1*k1, param=1)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1    0

Examples: Monomial Index-by-Order

julia> deriv(x1*k1, [1,0])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    1


julia> deriv(x1*k1, [0,1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1    0


julia> deriv(x1*k1, [1,1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    0        0    0

Examples: Monomial Index-by-Sparse Monomial

julia> deriv(x1*k1, [1=>1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    1


julia> deriv(x1*k1, params=[1=>1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1    0


julia> deriv(x1*k1, [1=>1], params=[1=>1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    0        0    0
GTPSA.evaluateMethod
evaluate(ct1::ComplexTPS, x::Vector{<:Number})::ComplexF64

Evaluates ct1 at the point x.

GTPSA.evaluateMethod
evaluate(t1::TPS, x::Vector{<:Real})::Float64

Evaluates t1 at the point x.

GTPSA.evaluateMethod
evaluate(m::Vector{ComplexTPS}, x::Vector{<:Number})::Vector{ComplexF64}

Evaluates m at the point x.

GTPSA.evaluateMethod
evaluate(m::Vector{TPS}, x::Vector{<:Real})::Vector{Float64}

Evaluates m at the point x.

GTPSA.exppbFunction
exppb(F::Vector{<:Union{TPS,ComplexTPS}}, m::Vector{<:Union{TPS,ComplexTPS}}=vars(first(F)))

Calculates exp(F⋅∇)m = m + F⋅∇m + (F⋅∇)²m/2! + .... If m is not provided, it is assumed to be the identity.

Example

julia> d = Descriptor(2,10); x = vars()[1]; p = vars()[2];

julia> time = 0.01; k = 2; m = 0.01;

julia> h = p^2/(2m) + 1/2*k*x^2;

julia> hf = getvectorfield(h);

julia> map = exppb(-time*hf, [x, p])
2-element Vector{TPS}:
  Out  Coefficient                Order   Exponent
-------------------------------------------------
   1:   9.9001665555952290e-01      1      1   0
   1:   9.9666999841313930e-01      1      0   1
-------------------------------------------------
   2:  -1.9933399968262787e-02      1      1   0
   2:   9.9001665555952378e-01      1      0   1
GTPSA.fgradMethod
fgrad(F::Vector{<:Union{TPS,ComplexTPS}}, g::Union{TPS,ComplexTPS})

Calculates F⋅∇g.

GTPSA.gethamiltonianMethod
gethamiltonian(F::Vector{<:Union{TPS,ComplexTPS}})

Assuming the variables in the TPSA are canonically-conjugate, and ordered so that the canonically- conjugate variables are consecutive (q₁, p₁, q₂, p₂, ...), this function calculates the Hamiltonian from a vector field F that can be obtained from a Hamiltonian (e.g. by getvectorfield). Explicitly, ∫ F₁ dp₁ - ∫ F₂ dq₁ + ... + ∫ F₂ₙ₋₁ dpₙ - ∫ F₂ₙ dqₙ

Example

julia> d = Descriptor(2,10); x = vars();

julia> h = (x[1]^2 + x[2]^2)/2
TPS:
 Coefficient                Order   Exponent
  5.0000000000000000e-01      2      2   0
  5.0000000000000000e-01      2      0   2


julia> F = getvectorfield(h)
2-element Vector{TPS}:
  Out  Coefficient                Order   Exponent
-------------------------------------------------
   1:  -1.0000000000000000e+00      1      0   1
-------------------------------------------------
   2:   1.0000000000000000e+00      1      1   0


julia> gethamiltonian(F)
TPS:
 Coefficient                Order   Exponent
  5.0000000000000000e-01      2      2   0
  5.0000000000000000e-01      2      0   2
GTPSA.getordMethod
getord(t1::Union{TPS, ComplexTPS}, order::Integer)

Extracts one homogenous polynomial from t1 of the given order.

GTPSA.getvectorfieldMethod
getvectorfield(h::Union{TPS,ComplexTPS})::Vector{<:typeof(h)}

Assuming the variables in the TPSA are canonically-conjugate, and ordered so that the canonically- conjugate variables are consecutive (q₁, p₁, q₂, p₂, ...), calculates the vector field (Hamilton's equations) from the passed Hamiltonian, defined as [∂h/∂p₁, -∂h/∂q₁, ...]

Example

julia> d = Descriptor(2,10); x = vars();

julia> h = (x[1]^2 + x[2]^2)/2
TPS:
 Coefficient                Order   Exponent
  5.0000000000000000e-01      2      2   0
  5.0000000000000000e-01      2      0   2


julia> getvectorfield(h)
2-element Vector{TPS}:
  Out  Coefficient                Order   Exponent
-------------------------------------------------
   1:  -1.0000000000000000e+00      1      0   1
-------------------------------------------------
   2:   1.0000000000000000e+00      1      1   0
GTPSA.gradient!Method
gradient!(result::Vector{ComplexF64}, ct::ComplexTPS; include_params=false)

Extracts the first-order partial derivatives (evaluated at 0) from the ComplexTPS and fills the result vector in-place. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the ComplexTPS.

Input

  • ctComplexTPS to extract the gradient from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • result – Preallocated Vector to fill with the gradient of the ComplexTPS
GTPSA.gradient!Method
gradient!(result::Vector{Float64}, t::TPS; include_params=false)

Extracts the first-order partial derivatives (evaluated at 0) from the TPS and fills the result vector in-place. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the TPS.

Input

  • tTPS to extract the gradient from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • result – Preallocated Vector to fill with the gradient of the TPS
GTPSA.gradientMethod
gradient(ct::ComplexTPS; include_params=false)::Vector{ComplexF64}

Extracts the first-order partial derivatives (evaluated at 0) from the ComplexTPS. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the ComplexTPS.

Input

  • ctComplexTPS to extract the gradient from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • grad – Gradient of the ComplexTPS
GTPSA.gradientMethod
gradient(t::TPS; include_params=false)::Vector{Float64}

Extracts the first-order partial derivatives (evaluated at 0) from the TPS. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the TPS.

Input

  • tTPS to extract the gradient from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • grad – Gradient of the TPS
GTPSA.hessian!Method
hessian!(result::Matrix{ComplexF64},ct::ComplexTPS; include_params=false)

Extracts the second-order partial derivatives (evaluated at 0) from the ComplexTPS and fills the result matrix in-place. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the second-order monomial coefficients already in the ComplexTPS.

Input

  • ctComplexTPS to extract the Hessian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • result – Preallocated matrix to fill with the Hessian of the ComplexTPS
GTPSA.hessian!Method
hessian!(result::Matrix{Float64},t::TPS; include_params=false)

Extracts the second-order partial derivatives (evaluated at 0) from the TPS and fills the result matrix in-place. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the second-order monomial coefficients already in the TPS.

Input

  • tTPS to extract the Hessian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • result – Preallocated matrix to fill with the Hessian of the TPS
GTPSA.hessianMethod
hessian(ct::ComplexTPS; include_params=false)::Matrix{ComplexF64}

Extracts the second-order partial derivatives (evaluated at 0) from the ComplexTPS. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the second-order monomial coefficients already in the ComplexTPS.

Input

  • tComplexTPS to extract the Hessian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • H – Hessian of the ComplexTPS
GTPSA.hessianMethod
hessian(t::TPS; include_params=false)::Matrix{Float64}

Extracts the second-order partial derivatives (evaluated at 0) from the TPS. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the second-order monomial coefficients already in the TPS.

Input

  • tTPS to extract the Hessian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • H – Hessian of the TPS
GTPSA.integFunction
integ(t1::Union{TPS, ComplexTPS}, var::Integer=1)
∫(t1::Union{TPS, ComplexTPS}, var::Integer=1)

Integrates t1 wrt the variable var. Integration wrt parameters is not allowed, and integration wrt higher order monomials is not currently supported.

GTPSA.jacobian!Method
jacobian!(result::Matrix{ComplexF64}, m::Vector{ComplexTPS}; include_params=false)

Extracts the first-order partial derivatives (evaluated at 0) from the Vector of ComplexTPSs. and fills the result matrix in-place. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the ComplexTPSs.

Input

  • mVector of ComplexTPSs. to extract the Jacobian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • result – Preallocated matrix to fill with the Jacobian of m
GTPSA.jacobian!Method
jacobian!(result::Matrix{Float64}, m::Vector{TPS}; include_params=false)

Extracts the first-order partial derivatives (evaluated at 0) from the Vector of TPSs. and fills the result matrix in-place. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the TPSs.

Input

  • mVector of TPSs. to extract the Jacobian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • result – Preallocated matrix to fill with the Jacobian of m
GTPSA.jacobianMethod
jacobian(m::Vector{ComplexTPS}; include_params=false)::Matrix{ComplexF64}

Extracts the first-order partial derivatives (evaluated at 0) from the Vector of ComplexTPSs. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the ComplexTPSs.

Input

  • mVector of ComplexTPSs. to extract the Jacobian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • J – Jacobian of m
GTPSA.jacobianMethod
jacobian(m::Vector{TPS}; include_params=false)::Matrix{Float64}

Extracts the first-order partial derivatives (evaluated at 0) from the Vector of TPSs. The partial derivatives wrt the parameters will also be extracted when the include_params flag is set to true. Note that this function is not calculating anything - just extracting the first-order monomial coefficients already in the TPSs.

Input

  • mVector of TPSs. to extract the Jacobian from
  • include_params – (Optional) Extract partial derivatives wrt parameters. Default is false

Output

  • J – Jacobian of m
GTPSA.lbMethod
lb(A::Vector{<:Union{TPS,ComplexTPS}}, F::Vector{<:Union{TPS,ComplexTPS}})

Computes the Lie bracket of the vector functions A and F, defined over N variables as Σᵢᴺ Aᵢ (∂F/∂xᵢ) - Fᵢ (∂A/∂xᵢ)

Example

julia> d = Descriptor(2,10); x = vars();

julia> A = [-x[2], x[1]]
2-element Vector{TPS}:
  Out  Coefficient                Order   Exponent
-------------------------------------------------
   1:  -1.0000000000000000e+00      1      0   1
-------------------------------------------------
   2:   1.0000000000000000e+00      1      1   0


julia> F = [-x[1]^2, 2*x[1]*x[2]]
2-element Vector{TPS}:
  Out  Coefficient                Order   Exponent
-------------------------------------------------
   1:  -1.0000000000000000e+00      2      2   0
-------------------------------------------------
   2:   2.0000000000000000e+00      2      1   1


julia> lb(A,F)
2-element Vector{TPS}:
  Out  Coefficient                Order   Exponent
-------------------------------------------------
   1:   4.0000000000000000e+00      2      1   1
-------------------------------------------------
   2:   3.0000000000000000e+00      2      2   0
   2:  -2.0000000000000000e+00      2      0   2
GTPSA.logpbFunction
logpb(mf::Vector{<:Union{TPS,ComplexTPS}}, mi::Vector{<:Union{TPS,ComplexTPS}}=vars(first(F)))

Given a final map mf and initial map mi, this function calculates the vector field F such that mf=exp(F⋅∇)mi. If mi is not provided, it is assumed to be the identity.

julia> d = Descriptor(2,10); x = vars()[1]; p = vars()[2];

julia> time = 0.01; k = 2; m = 0.01;

julia> h = p^2/(2m) + 1/2*k*x^2;

julia> hf = getvectorfield(h);

julia> map = exppb(-time*hf);

julia> logpb(map) == -time*hf
true
GTPSA.mad_ctpsa_acc!Method
mad_ctpsa_acc!(a::Ptr{CTPSA}, v::ComplexF64, c::Ptr{CTPSA})

Adds a*v to TPSA c. Aliasing OK.

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c += v*a
GTPSA.mad_ctpsa_acc_r!Method
mad_ctpsa_acc_r!(a::Ptr{CTPSA}, v_re::Cdouble, v_im::Cdouble, c::Ptr{CTPSA})

Adds a*v to TPSA c. Aliasing OK. Without complex-by-value arguments.

Input

  • a – Source TPSA a
  • v_re – Real part of scalar with double precision
  • v_im – Imaginary part of scalar with double precision

Output

  • c – Destination TPSA c += v*a
GTPSA.mad_ctpsa_acos!Method
mad_ctpsa_acos!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the acos of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = acos(a)
GTPSA.mad_ctpsa_acosh!Method
mad_ctpsa_acosh!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the acosh of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = acosh(a)
GTPSA.mad_ctpsa_acot!Method
mad_ctpsa_acot!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the acot of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = acot(a)
GTPSA.mad_ctpsa_acoth!Method
mad_ctpsa_acoth!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the acoth of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = acoth(a)
GTPSA.mad_ctpsa_add!Method
mad_ctpsa_add!(a::Ptr{CTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination TPSA c = a + b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a + b
GTPSA.mad_ctpsa_addt!Method
mad_ctpsa_addt!(a::Ptr{CTPSA}, b::Ptr{RTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a + b (internal real-to-complex conversion).

Input

  • a – Source CTPSA a
  • b – Source RTPSA b

Output

  • c – Destination CTPSA c = a + b
GTPSA.mad_ctpsa_asin!Method
mad_ctpsa_asin!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the asin of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = asin(a)
GTPSA.mad_ctpsa_asinc!Method
mad_ctpsa_asinc!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the asinc(a) = asin(a)/a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = asinc(a) = asin(a)/a
GTPSA.mad_ctpsa_asinh!Method
mad_ctpsa_asinh!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the asinh of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = asinh(a)
GTPSA.mad_ctpsa_asinhc!Method
mad_ctpsa_asinhc!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the asinhc of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = asinhc(a)
GTPSA.mad_ctpsa_atan!Method
mad_ctpsa_atan!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the atan of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = atan(a)
GTPSA.mad_ctpsa_atanh!Method
mad_ctpsa_atanh!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the atanh of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = atanh(a)
GTPSA.mad_ctpsa_ax2pby2pcz2!Method
mad_ctpsa_ax2pby2pcz2!(a::ComplexF64, x::Ptr{CTPSA}, b::ComplexF64, y::Ptr{CTPSA}, c::ComplexF64, z::Ptr{CTPSA}, r::Ptr{CTPSA})

r = a*x^2 + b*y^2 + c*z^2

Input

  • a – Scalar a
  • x – TPSA x
  • b – Scalar b
  • y – TPSA y
  • c – Scalar c
  • z – TPSA z

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_ax2pby2pcz2_r!Method
mad_ctpsa_ax2pby2pcz2_r!(a_re::Cdouble, a_im::Cdouble, x::Ptr{CTPSA}, b_re::Cdouble, b_im::Cdouble, y::Ptr{CTPSA}, c_re::Cdouble, c_im::Cdouble, z::Ptr{CTPSA}, r::Ptr{CTPSA})

r = a*x^2 + b*y^2 + c*z^2. Same as mad_ctpsa_ax2pby2pcz2 without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • x – TPSA x
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b
  • y – TPSA y
  • c_re – Real part of Scalar c
  • c_im – Imag part of Scalar c
  • z – TPSA z

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axpb!Method
mad_ctpsa_axpb!(a::ComplexF64, x::Ptr{CTPSA}, b::ComplexF64, r::Ptr{CTPSA})

r = a*x + b

Input

  • a – Scalar a
  • x – TPSA x
  • b – Scalar b

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axpb_r!Method
mad_ctpsa_axpb_r!(a_re::Cdouble, a_im::Cdouble, x::Ptr{CTPSA}, b_re::Cdouble, b_im::Cdouble, r::Ptr{CTPSA})

r = a*x + b. Same as mad_ctpsa_axpb without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • x – TPSA x
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axpbypc!Method
mad_ctpsa_axpbypc!(a::ComplexF64, x::Ptr{CTPSA}, b::ComplexF64, y::Ptr{CTPSA}, c::ComplexF64, r::Ptr{CTPSA})

r = a*x+b*y+c

Input

  • a – Scalar a
  • x – TPSA x
  • b – Scalar b
  • y – TPSA y
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axpbypc_r!Method
mad_ctpsa_axpbypc_r!(a_re::Cdouble, a_im::Cdouble, x::Ptr{CTPSA}, b_re::Cdouble, b_im::Cdouble, y::Ptr{CTPSA}, c_re::Cdouble, c_im::Cdouble, r::Ptr{CTPSA})

r = a*x + b*y + c. Same as mad_ctpsa_axpbypc without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • x – TPSA x
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b
  • y – TPSA y
  • c_re – Real part of Scalar c
  • c_im – Imag part of Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axpsqrtbpcx2!Method
mad_ctpsa_axpsqrtbpcx2!(x::Ptr{CTPSA}, a::ComplexF64, b::ComplexF64, c::ComplexF64, r::Ptr{CTPSA})

r = a*x + sqrt(b + c*x^2)

Input

  • x – TPSA x
  • a – Scalar a
  • b – Scalar b
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axpsqrtbpcx2_r!Method
mad_ctpsa_axpsqrtbpcx2_r!(x::Ptr{CTPSA}, a_re::Cdouble, a_im::Cdouble, b_re::Cdouble, b_im::Cdouble, c_re::Cdouble, c_im::Cdouble, r::Ptr{CTPSA})

r = a*x + sqrt(b + c*x^2). Same as mad_ctpsa_axpsqrtbpcx2 without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b
  • c_re – Real part of Scalar c
  • c_im – Imag part of Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axypb!Method
mad_ctpsa_axypb!(a::ComplexF64, x::Ptr{CTPSA}, y::Ptr{CTPSA}, b::ComplexF64, r::Ptr{CTPSA})

r = a*x*y + b

Input

  • a – Scalar a
  • x – TPSA x
  • y – TPSA y
  • b – Scalar b

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axypb_r!Method
mad_ctpsa_axypb_r!(a_re::Cdouble, a_im::Cdouble, x::Ptr{CTPSA}, y::Ptr{CTPSA}, b_re::Cdouble, b_im::Cdouble, r::Ptr{CTPSA})

r = a*x*y + b. Same as mad_ctpsa_axypb without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • x – TPSA x
  • y – TPSA y
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axypbvwpc!Method
mad_ctpsa_axypbvwpc!(a::ComplexF64, x::Ptr{CTPSA}, y::Ptr{CTPSA}, b::ComplexF64, v::Ptr{CTPSA}, w::Ptr{CTPSA}, c::ComplexF64, r::Ptr{CTPSA})

r = a*x*y + b*v*w + c

Input

  • a – Scalar a
  • x – TPSA x
  • y – TPSA y
  • b – Scalar b
  • v – TPSA v
  • w – TPSA w
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axypbvwpc_r!Method
mad_ctpsa_axypbvwpc_r!(a_re::Cdouble, a_im::Cdouble, x::Ptr{CTPSA}, y::Ptr{CTPSA}, b_re::Cdouble, b_im::Cdouble, v::Ptr{CTPSA}, w::Ptr{CTPSA}, c_re::Cdouble, c_im::Cdouble, r::Ptr{CTPSA})

r = a*x*y + b*v*w + c. Same as mad_ctpsa_axypbvwpc without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • x – TPSA x
  • y – TPSA y
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b
  • v – TPSA v
  • w – TPSA w
  • c_re – Real part of Scalar c
  • c_im – Imag part of Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axypbzpc!Method
mad_ctpsa_axypbzpc!(a::ComplexF64, x::Ptr{CTPSA}, y::Ptr{CTPSA}, b::ComplexF64, z::Ptr{CTPSA}, c::ComplexF64, r::Ptr{CTPSA})

r = a*x*y + b*z + c

Input

  • a – Scalar a
  • x – TPSA x
  • y – TPSA y
  • b – Scalar b
  • z – TPSA z
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_axypbzpc_r!Method
mad_ctpsa_axypbzpc_r!(a_re::Cdouble, a_im::Cdouble, x::Ptr{CTPSA}, y::Ptr{CTPSA}, b_re::Cdouble, b_im::Cdouble, z::Ptr{CTPSA}, c_re::Cdouble, c_im::Cdouble, r::Ptr{CTPSA})

r = a*x*y + b*z + c. Same as mad_ctpsa_axypbzpc without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • x – TPSA x
  • y – TPSA y
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b
  • z – TPSA z
  • c_re – Real part of Scalar c
  • c_im – Imag part of Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_cabs!Method
mad_ctpsa_cabs!(t::Ptr{CTPSA}, r::Ptr{RTPSA})

Sets the RTPSA r equal to the aboslute value of CTPSA t. Specifically, the result contains a TPSA with the abs of all coefficients.

Input

  • t – Source CTPSA

Output

  • r – Destination RTPSA with r = |t|
GTPSA.mad_ctpsa_carg!Method
mad_ctpsa_carg!(t::Ptr{CTPSA}, r::Ptr{RTPSA})

Sets the RTPSA r equal to the argument (phase) of CTPSA t

Input

  • t – Source CTPSA

Output

  • r – Destination RTPSA with r = carg(t)
GTPSA.mad_ctpsa_clear!Method
mad_ctpsa_clear!(t::Ptr{CTPSA})

Clears the TPSA (reset to 0)

Input

  • t – Complex TPSA
GTPSA.mad_ctpsa_compose!Method
mad_ctpsa_compose!(na::Cint, ma::Vector{Ptr{CTPSA}}, nb::Cint, mb::Vector{Ptr{CTPSA}}, mc::Vector{Ptr{CTPSA}})

Composes two maps.

Input

  • na – Number of TPSAs in Map ma
  • ma – Map ma
  • nb – Number of TPSAs in Map mb
  • mb – Map mb

Output

  • mc – Composition of maps ma and mb
GTPSA.mad_ctpsa_conj!Method
mad_ctpsa_conj(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Calculates the complex conjugate of of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = conj(a)
GTPSA.mad_ctpsa_convert!Method
mad_ctpsa_convert!(t::Ptr{CTPSA}, r::Ptr{CTPSA}, n::Cint, t2r_::Vector{Cint}, pb::Cint)

General function to convert TPSAs to different orders and reshuffle canonical coordinates. The destination TPSA will be of order n, and optionally have the variable reshuffling defined by t2r_ and poisson bracket sign. e.g. if t2r_ = {1,2,3,4,6,5} and pb = -1, canonical coordinates 6 and 5 are swapped and the new 5th canonical coordinate will be negated. Useful for comparing with different differential algebra packages.

Input

  • t – Source complex TPSA
  • n – Length of vector
  • t2r_ – (Optional) Vector of index lookup
  • pb – Poisson bracket, 0, 1:fwd, -1:bwd

Output

  • r – Destination complex TPSA with specified order and canonical coordinate reshuffling.
GTPSA.mad_ctpsa_copy!Method
mad_ctpsa_copy!(t::Ptr{CTPSA}, r::Ptr{CTPSA})

Makes a copy of the complex TPSA t to r.

Input

  • t – Source complex TPSA

Output

  • r – Destination complex TPSA
GTPSA.mad_ctpsa_cos!Method
mad_ctpsa_cos!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the cos of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = cos(a)
GTPSA.mad_ctpsa_cosh!Method
mad_ctpsa_cosh!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the cosh of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = cosh(a)
GTPSA.mad_ctpsa_cot!Method
mad_ctpsa_cot!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the cot of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = cot(a)
GTPSA.mad_ctpsa_coth!Method
mad_ctpsa_coth!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the coth of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = coth(a)
GTPSA.mad_ctpsa_cplx!Method
mad_ctpsa_cplx!(re_::Ptr{RTPSA}, im_::Ptr{RTPSA}, r::Ptr{CTPSA})

Creates a CTPSA with real and imaginary parts from the RTPSAs re_ and im_ respectively.

Input

  • re_ – Real part of CTPSA to make
  • im_ – Imaginary part of CTPSA to make

Output

  • r – Destination CTPSA with r = re_ + im*im_
GTPSA.mad_ctpsa_cutord!Method
mad_ctpsa_cutord!(t::Ptr{CTPSA}, r::Ptr{CTPSA}, ord::Cint)

Cuts the TPSA off at the given order and above, or if ord is negative, will cut orders below abs(ord) (e.g. if ord = -3, then orders 0-3 are cut off).

Input

  • t – Source complex TPSA
  • ord – Cut order: 0..-ord or ord..mo

Output

  • r – Destination complex TPSA
GTPSA.mad_ctpsa_cycle!Method
mad_ctpsa_cycle!(t::Ptr{CTPSA}, i::Cint, n::Cint, m_::Vector{Cuchar}, v_::Ref{ComplexF64})::Cint

Used for scanning through each nonzero monomial in the TPSA. Given a starting index (-1 if starting at 0), will optionally fill monomial m_ with the monomial at index i and the value at v_ with the monomials coefficient, and return the next NONZERO monomial index in the TPSA. This is useful for building an iterator through the TPSA.

Input

  • t – TPSA to scan
  • i – Index to start from (-1 to start at 0)
  • n – Size of monomial
  • m_ – (Optional) Monomial to be filled if provided
  • v_ – (Optional) Pointer to value of coefficient

Output

  • i – Index of next nonzero monomial in the TPSA, or -1 if reached the end
GTPSA.mad_ctpsa_debugMethod
mad_ctpsa_debug(t::Ptr{CTPSA}, name_::Cstring, fnam_::Cstring, line_::Cint, stream_::Ptr{Cvoid})

Prints TPSA with all information of data structure.

Input

  • t – TPSA
  • name_ – (Optional) Name of TPSA
  • fnam_ – (Optional) File name to print to
  • line_ – (Optional) Line number in file to start at
  • stream_ – (Optional) I/O stream to print to, default is stdout
GTPSA.mad_ctpsa_del!Method
mad_ctpsa_del!(t::Ptr{CTPSA})

Calls the destructor for the complex TPSA.

Input

  • t – Complex TPSA to destruct
GTPSA.mad_ctpsa_deriv!Method
mad_ctpsa_deriv!(a::Ptr{CTPSA}, c::Ptr{CTPSA}, iv::Cint)

Differentiates TPSA with respect to the variable with index iv.

Input

  • a – Source TPSA to differentiate
  • iv – Index of variable to take derivative wrt to (e.g. derivative wrt x, iv = 1).

Output

  • c – Destination TPSA
GTPSA.mad_ctpsa_derivm!Method
mad_ctpsa_derivm!(a::Ptr{CTPSA}, c::Ptr{CTPSA}, n::Cint, m::Vector{Cuchar})

Differentiates TPSA with respect to the monomial defined by byte array m.

Input

  • a – Source TPSA to differentiate
  • n – Length of monomial to differentiate wrt
  • m – Monomial to take derivative wrt

Output

  • c – Destination TPSA
GTPSA.mad_ctpsa_descMethod
mad_ctpsa_desc(t::Ptr{CTPSA})::Ptr{Desc}

Gets the descriptor for the complex TPSA.

Input

  • t – Complex TPSA

Output

  • ret – Descriptor for the TPSA
GTPSA.mad_ctpsa_dif!Method
mad_ctpsa_dif!(a::Ptr{CTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

For each homogeneous polynomial in TPSAs a and b, calculates either the relative error or absolute error for each order. If the maximum coefficient for a given order in a is > 1, the relative error is computed for that order. Else, the absolute error is computed. This is very useful for comparing maps between codes or doing unit tests. In Julia, essentially:

c_i = (a_i.-b_i)/maximum([abs.(a_i)...,1]) where a_i and b_i are vectors of the monomials for an order i

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c
GTPSA.mad_ctpsa_dift!Method
mad_ctpsa_dift!(a::Ptr{CTPSA}, b::Ptr{RTPSA}, c::Ptr{CTPSA})

For each homogeneous polynomial in CTPSA a and RTPSA b, calculates either the relative error or absolute error for each order. If the maximum coefficient for a given order in a is > 1, the relative error is computed for that order. Else, the absolute error is computed. This is very useful for comparing maps between codes or doing unit tests. In Julia, essentially:

c_i = (a_i.-b_i)/maximum([abs.(a_i)...,1]) where a_i and b_i are vectors of the monomials for an order i

Input

  • a – Source CTPSA a
  • b – Source RTPSA b

Output

  • c – Destination CTPSA c
GTPSA.mad_ctpsa_div!Method
mad_ctpsa_div!(a::Ptr{CTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination TPSA c = a / b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a / b
GTPSA.mad_ctpsa_divt!Method
mad_ctpsa_divt!(a::Ptr{CTPSA}, b::Ptr{RTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a / b (internal real-to-complex conversion).

Input

  • a – Source CTPSA a
  • b – Source RTPSA b

Output

  • c – Destination CTPSA c = a / b
GTPSA.mad_ctpsa_equMethod
mad_ctpsa_equ(a::Ptr{CTPSA}, b::Ptr{CTPSA}, tol_::Cdouble)::Cuchar

Checks if the TPSAs a and b are equal within the specified tolerance tol_. If tol_ is not specified, DBL_GTPSA.show_epsILON is used.

Input

  • a – TPSA a
  • b – TPSA b
  • tol_ – (Optional) Difference below which the TPSAs are considered equal

Output

  • ret - True if a == b within tol_
GTPSA.mad_ctpsa_equtMethod
mad_ctpsa_equt(a::Ptr{CTPSA}, b::Ptr{RTPSA}, tol_::Cdouble)::Cuchar

Checks if the CTPSA a is equal to the RTPSA b within the specified tolerance tol_ (internal real-to-complex conversion).

Input

  • a – CTPSA a
  • b – RTPSA b
  • tol_ – (Optional) Difference below which the TPSAs are considered equal

Output

  • ret - True if a == b within tol_
GTPSA.mad_ctpsa_erf!Method
mad_ctpsa_erf!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the erf of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = erf(a)
GTPSA.mad_ctpsa_erfc!Method
mad_ctpsa_erfc!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the erfc of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = erfc(a)
GTPSA.mad_ctpsa_eval!Method
mad_ctpsa_eval!(na::Cint, ma::Vector{Ptr{CTPSA}}, nb::Cint, tb::Vector{ComplexF64}, tc::Vector{ComplexF64})

Evaluates the map at the point tb

Input

  • na – Number of TPSAs in the map
  • ma – Map ma
  • nb – Length of tb
  • tb – Point at which to evaluate the map

Output

  • tc – Values for each TPSA in the map evaluated at the point tb
GTPSA.mad_ctpsa_exp!Method
mad_ctpsa_exp!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the exp of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = exp(a)
GTPSA.mad_ctpsa_exppb!Method
mad_ctpsa_exppb!(na::Cint, ma::Vector{Ptr{CTPSA}}, mb::Vector{Ptr{CTPSA}}, mc::Vector{Ptr{CTPSA}})

Computes the exponential of fgrad of the vector fields ma and mb, literally exppb(ma, mb) = mb + fgrad(ma, mb) + fgrad(ma, fgrad(ma, mb))/2! + ...

Input

  • na – Length of ma and mb
  • ma – Vector of TPSA ma
  • mb – Vector of TPSA mb

Output

  • mc – Destination vector of TPSA mc
GTPSA.mad_ctpsa_fgrad!Method
mad_ctpsa_fgrad!(na::Cint, ma::Vector{Ptr{CTPSA}}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Calculates dot(ma, grad(b))

Input

  • na – Length of ma consistent with number of variables in b
  • ma – Vector of TPSA
  • b – TPSA

Output

  • cdot(ma, grad(b))
GTPSA.mad_ctpsa_fld2vec!Method
mad_ctpsa_fld2vec!(na::Cint, ma::Vector{Ptr{CTPSA}}, c::Ptr{CTPSA})

Assuming the variables in the TPSA are canonically-conjugate, and ordered so that the canonically- conjugate variables are consecutive (q1, p1, q2, p2, ...), calculates the Hamiltonian one obtains from ther vector field (in the form [da/dp1, -da/dq1, ...])

Input

  • na – Number of TPSA in ma consistent with number of variables in c
  • ma – Vector field

Output

  • c – Hamiltonian as a TPSA derived from the vector field ma
GTPSA.mad_ctpsa_get0Method
mad_ctpsa_get0(t::Ptr{CTPSA})::ComplexF64

Gets the 0th order (scalar) value of the TPSA

Input

  • t – TPSA

Output

  • ret – Scalar value of TPSA
GTPSA.mad_ctpsa_get0_r!Method
mad_ctpsa_get0_r!(t::Ptr{CTPSA}, r::Ref{ComplexF64})

Gets the 0th order (scalar) value of the TPSA in place.

Input

  • t – TPSA

Output

  • r – Scalar value of TPSA
GTPSA.mad_ctpsa_getiMethod
mad_ctpsa_geti(t::Ptr{CTPSA}, i::Cint)::ComplexF64

Gets the coefficient of the monomial at index i. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • i – Monomial index

Output

  • ret – Coefficient of monomial at index i
GTPSA.mad_ctpsa_geti_r!Method
mad_ctpsa_geti_r!(t::Ptr{CTPSA}, i::Cint,  r::Ref{ComplexF64})

Gets the coefficient of the monomial at index i in place. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • i – Monomial index

Output

  • r – Coefficient of monomial at index i
GTPSA.mad_ctpsa_getmMethod
mad_ctpsa_getm(t::Ptr{CTPSA}, n::Cint, m::Vector{Cuchar})::ComplexF64

Gets the coefficient of the monomial m defined as a byte array. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as byte array

Output

  • ret – Coefficient of monomial m in TPSA
GTPSA.mad_ctpsa_getm_r!Method
mad_ctpsa_getm_r!(t::Ptr{CTPSA}, n::Cint, m::Vector{Cuchar}, r::Ref{ComplexF64})

Gets the coefficient of the monomial m defined as a byte array in place. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as byte array

Output

  • r – Coefficient of monomial m in TPSA
GTPSA.mad_ctpsa_getord!Method
mad_ctpsa_getord!(t::Ptr{CTPSA}, r::Ptr{CTPSA}, ord::Cuchar)

Extract one homogeneous polynomial of the given order

Input

  • t – Sourcecomplex TPSA
  • ord – Order to retrieve

Output

  • r – Destination complex TPSA
GTPSA.mad_ctpsa_getsMethod
mad_ctpsa_gets(t::Ptr{CTPSA}, n::Cint, s::Cstring)::ComplexF64

Gets the coefficient of the monomial s defined as a string. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Size of monomial
  • s – Monomial as string

Output

  • ret – Coefficient of monomial s in TPSA
GTPSA.mad_ctpsa_gets_r!Method
mad_ctpsa_gets_r!(t::Ptr{CTPSA}, n::Cint, s::Cstring, r::Ref{ComplexF64})

Gets the coefficient of the monomial s defined as a string in place. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as string

Output

  • r – Coefficient of monomial s in TPSA
GTPSA.mad_ctpsa_getsmMethod
mad_ctpsa_getsm(t::Ptr{CTPSA}, n::Cint, m::Vector{Cint})::ComplexF64

Gets the coefficient of the monomial m defined as a sparse monomial. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as sparse monomial

Output

  • ret – Coefficient of monomial m in TPSA
GTPSA.mad_ctpsa_getsm_r!Method
mad_ctpsa_getsm_r!(t::Ptr{CTPSA}, n::Cint, m::Vector{Cint}, r::Ref{ComplexF64})

Gets the coefficient of the monomial m defined as a sparse monomial in place. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as sparse monomial

Output

  • r – Coefficient of monomial m in TPSA
GTPSA.mad_ctpsa_getv!Method
mad_ctpsa_getv!(t::Ptr{CTPSA}, i::Cint, n::Cint, v::Vector{ComplexF64})

Vectorized getter of the coefficients for monomials with indices i..i+n. Useful for extracting the 1st order parts of a TPSA to construct a matrix (i = 1, n = nv+np = nn).

Input

  • t – TPSA
  • i – Starting index of monomials to get coefficients
  • n – Number of monomials to get coefficients of starting at i

Output

  • v – Array of coefficients for monomials i..i+n
GTPSA.mad_ctpsa_hypot!Method
mad_ctpsa_hypot!(x::Ptr{CTPSA}, y::Ptr{CTPSA}, r::Ptr{CTPSA})

Sets TPSA r to sqrt(real(x)^2+real(y)^2) + im*sqrt(imag(x)^2+imag(y)^2)

Input

  • x – Source TPSA x
  • y – Source TPSA y

Output

  • r – Destination TPSA sqrt(real(x)^2+real(y)^2) + im*sqrt(imag(x)^2+imag(y)^2)
GTPSA.mad_ctpsa_hypot3!Method
mad_ctpsa_hypot3!(x::Ptr{CTPSA}, y::Ptr{CTPSA}, z::Ptr{CTPSA}, r::Ptr{CTPSA})

Sets TPSA r to sqrt(x^2+y^2+z^2). Does NOT allow for r = x, y, z !!!

Input

  • x – Source TPSA x
  • y – Source TPSA y
  • z – Source TPSA z

Output

  • r – Destination TPSA r = sqrt(x^2+y^2+z^2)
GTPSA.mad_ctpsa_idxmMethod
mad_ctpsa_idxm(t::Ptr{CTPSA}, n::Cint, m::Vector{Cuchar})::Cint

Returns index of monomial in the TPSA given the monomial as a byte array. This generally should not be used, as there are no assumptions about which monomial is attached to which index.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as byte array

Output

  • ret – Index of monomial in TPSA
GTPSA.mad_ctpsa_idxsMethod
mad_ctpsa_idxs(t::Ptr{CTPSA}, n::Cint, s::Cstring)::Cint

Returns index of monomial in the TPSA given the monomial as string. This generally should not be used, as there are no assumptions about which monomial is attached to which index.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as string

Output

  • ret – Index of monomial in TPSA
GTPSA.mad_ctpsa_idxsmMethod
mad_ctpsa_idxsm(t::Ptr{CTPSA}, n::Cint, m::Vector{Cint})::Cint

Returns index of monomial in the TPSA given the monomial as a sparse monomial. This generally should not be used, as there are no assumptions about which monomial is attached to which index.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as sparse monomial

Output

  • ret – Index of monomial in TPSA
GTPSA.mad_ctpsa_imag!Method
mad_ctpsa_imag!(t::Ptr{CTPSA}, r::Ptr{RTPSA})

Sets the RTPSA r equal to the imaginary part of CTPSA t.

Input

  • t – Source CTPSA

Output

  • r – Destination RTPSA with r = Im(t)
GTPSA.mad_ctpsa_init!Method
mad_ctpsa_init(t::Ptr{CTPSA}, d::Ptr{Desc}, mo::Cuchar)::Ptr{CTPSA}

Unsafe initialization of an already existing TPSA t with maximum order mo to the descriptor d. mo must be less than the maximum order of the descriptor. t is modified in place and also returned.

Input

  • t – TPSA to initialize to descriptor d
  • d – Descriptor
  • mo – Maximum order of the TPSA (must be less than maximum order of the descriptor)

Output

  • t – TPSA initialized to descriptor d with maximum order mo
GTPSA.mad_ctpsa_integ!Method
mad_ctpsa_integ!(a::Ptr{CTPSA}, c::Ptr{CTPSA}, iv::Cint)

Integrates TPSA with respect to the variable with index iv.

Input

  • a – Source TPSA to integrate
  • iv – Index of variable to integrate over (e.g. integrate over x, iv = 1).

Output

  • c – Destination TPSA
GTPSA.mad_ctpsa_inv!Method
mad_ctpsa_inv!(a::Ptr{CTPSA},  v::ComplexF64, c::Ptr{CTPSA})

Sets TPSA c to v/a.

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c = v/a
GTPSA.mad_ctpsa_inv_r!Method
mad_ctpsa_inv_r!(a::Ptr{CTPSA}, v_re::Cdouble, v_im::Cdouble, c::Ptr{CTPSA})

Sets TPSA c to v/a. Without complex-by-value arguments.

Input

  • a – Source TPSA a
  • v_re – Real part of scalar with double precision
  • v_im – Imaginary part of scalar with double precision

Output

  • c – Destination TPSA c = v*a
GTPSA.mad_ctpsa_invsqrt!Method
mad_ctpsa_invsqrt!(a::Ptr{CTPSA}, v::ComplexF64, c::Ptr{CTPSA})

Sets TPSA c to v/sqrt(a).

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c = v/sqrt(a)
GTPSA.mad_ctpsa_invsqrt_r!Method
mad_ctpsa_invsqrt_r!(a::Ptr{CTPSA}, v_re::Cdouble, v_im::Cdouble, c::Ptr{CTPSA})

Sets TPSA c to v/sqrt(a). Without complex-by-value arguments.

Input

  • a – Source TPSA a
  • v_re – Real part of scalar with double precision
  • v_im – Imaginary part of scalar with double precision

Output

  • c – Destination TPSA c = v*a
GTPSA.mad_ctpsa_isnulMethod
mad_ctpsa_isnul(t::Ptr{CTPSA})::Cuchar

Checks if TPSA is 0 or not

Input

  • t – Complex TPSA to check

Output

  • ret – True or false
GTPSA.mad_ctpsa_isvalidMethod
mad_ctpsa_isvalid(t::Ptr{CTPSA})::Cuchar

Sanity check of the TPSA integrity.

Input

  • t – Complex TPSA to check if valid

Output

  • ret – True if valid TPSA, false otherwise
GTPSA.mad_ctpsa_lenMethod
mad_ctpsa_len(t::Ptr{CTPSA})::Cint

Gets the length of the TPSA itself (e.g. the descriptor may be order 10 but TPSA may only be order 2)

Input

  • t – Complex TPSA

Output

  • ret – Length of CTPSA
GTPSA.mad_ctpsa_liebra!Method
mad_ctpsa_liebra!(na::Cint, ma::Vector{Ptr{CTPSA}}, mb::Vector{Ptr{CTPSA}}, mc::Vector{Ptr{CTPSA}})

Computes the Lie bracket of the vector fields ma and mb, defined as sumi mai (dmb/dxi) - mbi (dma/dx_i).

Input

  • na – Length of ma and mb
  • ma – Vector of TPSA ma
  • mb – Vector of TPSA mb

Output

  • mc – Destination vector of TPSA mc
GTPSA.mad_ctpsa_log!Method
mad_ctpsa_log!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the log of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = log(a)
GTPSA.mad_ctpsa_logaxpsqrtbpcx2!Method
mad_ctpsa_logaxpsqrtbpcx2!(x::Ptr{CTPSA}, a::ComplexF64, b::ComplexF64, c::ComplexF64, r::Ptr{CTPSA})

r = log(a*x + sqrt(b + c*x^2))

Input

  • x – TPSA x
  • a – Scalar a
  • b – Scalar b
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_logaxpsqrtbpcx2_r!Method
mad_ctpsa_logaxpsqrtbpcx2_r!(x::Ptr{CTPSA}, a_re::Cdouble, a_im::Cdouble, b_re::Cdouble, b_im::Cdouble, c_re::Cdouble, c_im::Cdouble, r::Ptr{CTPSA})

r = log(a*x + sqrt(b + c*x^2)). Same as mad_ctpsa_logaxpsqrtbpcx2 without complex-by-value arguments.

Input

  • a_re – Real part of Scalar a
  • a_im – Imag part of Scalar a
  • b_re – Real part of Scalar b
  • b_im – Imag part of Scalar b
  • c_re – Real part of Scalar c
  • c_im – Imag part of Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_logpb!Method
mad_ctpsa_logpb!(na::Cint, ma::Vector{Ptr{CTPSA}}, mb::Vector{Ptr{CTPSA}}, mc::Vector{Ptr{CTPSA}})

Computes the log of the Poisson bracket of the vector of TPSA ma and mb; the result is the vector field F used to evolve to ma from mb.

Input

  • na – Length of ma and mb
  • ma – Vector of TPSA ma
  • mb – Vector of TPSA mb

Output

  • mc – Destination vector of TPSA mc
GTPSA.mad_ctpsa_logxdy!Method
mad_ctpsa_logxdy!(x::Ptr{CTPSA}, y::Ptr{CTPSA}, r::Ptr{CTPSA})

r = log(x / y)

Input

  • x – TPSA x
  • y – TPSA y

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_maxordMethod
mad_ctpsa_maxord(t::Ptr{CTPSA}, n::Cint, idx_::Vector{Cint})::Cint

Returns the index to the monomial with maximum abs(coefficient) in the TPSA for all orders 0 to n. If idx_ is provided, it is filled with the indices for the maximum abs(coefficient) monomial for each order up to n.

Input

  • t – Complex TPSA
  • n – Highest order to include in finding the maximum abs(coefficient) in the TPSA, length of idx_ if provided

Output

  • idx_ – (Optional) If provided, is filled with indices to the monomial for each order up to n with maximum abs(coefficient)
  • mi – Index to the monomial in the TPSA with maximum abs(coefficient)
GTPSA.mad_ctpsa_mconv!Method
mad_ctpsa_mconv!(na::Cint, ma::Vector{Ptr{CTPSA}}, nc::Cint, mc::Vector{Ptr{CTPSA}}, n::Cint, t2r_::Vector{Cint}, pb::Cint)

Equivalent to mad_tpsa_convert, but applies the conversion to all TPSAs in the map ma.

Input

  • na – Number of TPSAs in the map
  • ma – Map ma
  • nc – Number of TPSAs in the output map mc
  • n – Length of vector (size of t2r_)
  • t2r_ – (Optional) Vector of index lookup
  • pb – Poisson bracket, 0, 1:fwd, -1:bwd

Output

  • mc – Map mc with specified conversions
GTPSA.mad_ctpsa_minv!Method
mad_ctpsa_minv!(na::Cint, ma::Vector{Ptr{CTPSA}}, mc::Vector{Ptr{CTPSA}})

Inverts the map.

Input

  • na – Number of TPSAs in the map
  • ma – Map ma

Output

  • mc – Inversion of Map ma
GTPSA.mad_ctpsa_mnrmMethod
mad_ctpsa_mnrm(na::Cint, ma::Vector{Ptr{CTPSA}})::Cdouble

Computes the norm of the map (sum of absolute value of coefficients of all TPSAs in the map).

Input

  • na – Number of TPSAs in the map
  • ma – Map ma

Output

  • nrm – Norm of map (sum of absolute value of coefficients of all TPSAs in the map)
GTPSA.mad_ctpsa_mono!Method
mad_ctpsa_mono!(t::Ptr{CTPSA}, i::Cint, n::Cint, m_::Vector{Cuchar}, p_::Vector{Cuchar})::Cuchar

Returns the order of the monomial at index i in the TPSA and optionally the monomial at that index is returned in m_ and the order of parameters in the monomial in p_

Input

  • t – TPSA
  • i – Index valid in TPSA
  • n – Length of monomial

Output

  • m_ – (Optional) Monomial at index i in TPSA
  • p_ – (Optional) Order of parameters in monomial
  • ret – Order of monomial in TPSA a index i
GTPSA.mad_ctpsa_mul!Method
mad_ctpsa_mul!(a::Ptr{CTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination TPSA c = a * b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a * b
GTPSA.mad_ctpsa_mult!Method
mad_ctpsa_mult!(a::Ptr{CTPSA}, b::Ptr{RTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a * b (internal real-to-complex conversion).

Input

  • a – Source CTPSA a
  • b – Source RTPSA b

Output

  • c – Destination CTPSA c = a * b
GTPSA.mad_ctpsa_namMethod
mad_ctpsa_nam(t::Ptr{CTPSA})::Cstring

Get the name of the TPSA.

Input

  • t – Complex TPSA

Output

  • ret – Name of CTPSA (Null terminated in C)
GTPSA.mad_ctpsa_newMethod
mad_ctpsa_new(t::Ptr{CTPSA}, mo::Cuchar)::Ptr{CTPSA}

Creates a blank TPSA with same number of variables/parameters of the inputted TPSA, with maximum order specified by mo. If MAD_TPSA_SAME is passed for mo, the mo currently in t is used for the created TPSA. Ok with t=(tpsa_t*)ctpsa

Input

  • t – TPSA
  • mo – Maximum order of new TPSA

Output

  • ret – New blank TPSA with maximum order mo
GTPSA.mad_ctpsa_newdMethod
mad_ctpsa_newd(d::Ptr{Desc}, mo::Cuchar)::Ptr{CTPSA}

Creates a complex TPSA defined by the specified descriptor and maximum order. If MADCTPSADEFAULT is passed for mo, the mo defined in the descriptor is used. If mo > d_mo, then mo = d_mo.

Input

  • d – Descriptor
  • mo – Maximum order

Output

  • t – New complex TPSA defined by the descriptor
GTPSA.mad_ctpsa_nrmMethod
mad_ctpsa_nrm(a::Ptr{CTPSA})::Cdouble

Calculates the 1-norm of TPSA a (sum of abs of all coefficients)

Input

  • a – TPSA

Output

  • nrm – 1-Norm of TPSA a
GTPSA.mad_ctpsa_ordMethod
mad_ctpsa_ord(t::Ptr{CTPSA})::Cuchar

Gets the TPSA order.

Input

  • t – Complex TPSA

Output

  • ret – Order of TPSA
GTPSA.mad_ctpsa_ordnMethod
mad_ctpsa_ordn(n::Cint, t::Vector{Ptr{CTPSA}})::Cuchar

Returns the max order of all TPSAs in t.

Input

  • n – Number of TPSAs
  • t – Array of TPSAs

Output

  • mo – Maximum order of all TPSAs
GTPSA.mad_ctpsa_ordvMethod
mad_ctpsa_ordv(t::Ptr{CTPSA}, ts::Ptr{CTPSA}...)::Cuchar

Returns maximum order of all TPSAs provided.

Input

  • t – TPSA
  • ts – Variable number of TPSAs passed as parameters

Output

  • mo – Maximum order of all TPSAs provided
GTPSA.mad_ctpsa_pminv!Method
mad_ctpsa_pminv!(na::Cint, ma::Vector{Ptr{CTPSA}}, mc::Vector{Ptr{CTPSA}}, select::Vector{Cint})

Computes the partial inverse of the map with only the selected variables, specified by 0s or 1s in select.

Input

  • na – Number of TPSAs in ma
  • ma – Map ma
  • select – Array of 0s or 1s defining which variables to do inverse on (atleast same size as na)

Output

  • mc – Partially inverted map using variables specified as 1 in the select array
GTPSA.mad_ctpsa_poisbra!Method
mad_ctpsa_poisbra!(a::Ptr{CTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA}, nv::Cint)

Sets TPSA c to the poisson bracket of TPSAs a and b.

Input

  • a – Source TPSA a
  • b – Source TPSA b
  • nv – Number of variables in the TPSA

Output

  • c – Destination TPSA c
GTPSA.mad_ctpsa_poisbrat!Method
mad_ctpsa_poisbrat!(a::Ptr{CTPSA}, b::Ptr{RTPSA}, c::Ptr{CTPSA}, nv::Cint)

Sets TPSA c to the poisson bracket of CTPSA aand RTPSA b (internal real-to-complex conversion).

Input

  • a – Source CTPSA a
  • b – Source RTPSA b
  • nv – Number of variables in the TPSA

Output

  • c – Destination CTPSA c
GTPSA.mad_ctpsa_polar!Method
mad_ctpsa_polar!(t::Ptr{CTPSA}, r::Ptr{CTPSA})

Sets r = |t| + im*atan2(Im(t), Re(t))

Input

  • t – Source CTPSA
  • r – Destination CTPSA
GTPSA.mad_ctpsa_pow!Method
mad_ctpsa_pow!(a::Ptr{CTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination TPSA c = a ^ b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a ^ b
GTPSA.mad_ctpsa_powi!Method
mad_ctpsa_powi!(a::Ptr{CTPSA}, n::Cint, c::Ptr{CTPSA})

Sets the destination TPSA c = a ^ n where n is an integer.

Input

  • a – Source TPSA a
  • n – Integer power

Output

  • c – Destination TPSA c = a ^ n
GTPSA.mad_ctpsa_pown!Method
mad_ctpsa_pown!(a::Ptr{CTPSA}, v::ComplexF64, c::Ptr{CTPSA})

Sets the destination TPSA c = a ^ v where v is of double precision.

Input

  • a – Source TPSA a
  • v – Power, ComplexF64

Output

  • c – Destination TPSA c = a ^ v
GTPSA.mad_ctpsa_pown_r!Method
mad_ctpsa_pown_r!(a::Ptr{CTPSA}, v_re::Cdouble, v_im::Cdouble, c::Ptr{CTPSA})

Sets the destination TPSA c = a ^ v where v is of double precision. Without complex-by-value arguments.

Input

  • a – Source TPSA a
  • v_re – Real part of power
  • v_im – Imaginary part of power

Output

  • c – Destination TPSA c = a ^ v
GTPSA.mad_ctpsa_powt!Method
mad_ctpsa_powt!(a::Ptr{CTPSA}, b::Ptr{RTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a ^ b (internal real-to-complex conversion).

Input

  • a – Source CTPSA a
  • b – Source RTPSA b

Output

  • c – Destination CTPSA c = a ^ b
GTPSA.mad_ctpsa_printMethod
mad_ctpsa_print(t::Ptr{CTPSA}, name_::Cstring, eps_::Cdouble, nohdr_::Cint, stream_::Ptr{Cvoid})

Prints the TPSA coefficients with precision eps_. If nohdr_ is not zero, the header is not printed.

Input

  • t – TPSA to print
  • name_ – (Optional) Name of TPSA
  • eps_ – (Optional) Precision to output
  • nohdr_ – (Optional) If True, no header is printed
  • stream_ – (Optional) FILE pointer of output stream. Default is stdout
GTPSA.mad_ctpsa_real!Method
mad_ctpsa_real!(t::Ptr{CTPSA}, r::Ptr{RTPSA})

Sets the RTPSA r equal to the real part of CTPSA t.

Input

  • t – Source CTPSA

Output

  • r – Destination RTPSA with r = Re(t)
GTPSA.mad_ctpsa_rect!Method
mad_ctpsa_rect!(t::Ptr{CTPSA}, r::Ptr{CTPSA})

Sets r = Re(t)*cos(Im(t)) + im*Re(t)*sin(Im(t))

Input

  • t – Source CTPSA
  • r – Destination CTPSA
GTPSA.mad_ctpsa_scanMethod
mad_ctpsa_scan(stream_::Ptr{Cvoid})::Ptr{CTPSA}

Scans in a TPSA from the stream_.

Input

  • stream_ – (Optional) I/O stream from which to read the TPSA, default is stdin

Output

  • t – TPSA scanned from I/O stream_
GTPSA.mad_ctpsa_scan_coef!Method
mad_ctpsa_scan_coef!(t::Ptr{CTPSA}, stream_::Ptr{Cvoid})

Read TPSA coefficients into TPSA t. This should be used with mad_tpsa_scan_hdr for external languages using this library where the memory is managed NOT on the C side.

Input

  • stream_ – (Optional) I/O stream to read TPSA from, default is stdin

Output

  • t – TPSA with coefficients scanned from stream_
GTPSA.mad_ctpsa_scan_hdrMethod
mad_ctpsa_scan_hdr(kind_::Ref{Cint}, name_::Ptr{Cuchar}, stream_::Ptr{Cvoid})::Ptr{Desc}

Read TPSA header. Returns descriptor for TPSA given the header. This is useful for external languages using this library where the memory is managed NOT on the C side.

Input

  • kind_ – (Optional) Real or complex TPSA, or detect automatically if not provided.
  • name_ – (Optional) Name of TPSA
  • stream_ – (Optional) I/O stream to read TPSA from, default is stdin

Output

  • ret – Descriptor for the TPSA
GTPSA.mad_ctpsa_scl!Method
mad_ctpsa_scl!(a::Ptr{CTPSA}, v::ComplexF64, c::Ptr{CTPSA})

Sets TPSA c to v*a.

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c = v*a
GTPSA.mad_ctpsa_scl_r!Method
mad_ctpsa_scl_r!(a::Ptr{CTPSA}, v_re::Cdouble, v_im::Cdouble,, c::Ptr{CTPSA})

Sets TPSA c to v*a. Without complex-by-value arguments.

Input

  • a – Source TPSA a
  • v_re – Real part of scalar with double precision
  • v_im – Imaginary part of scalar with double precision

Output

  • c – Destination TPSA c = v*a
GTPSA.mad_ctpsa_sclord!Method
mad_ctpsa_sclord!(t::Ptr{CTPSA}, r::Ptr{CTPSA}, inv::Cuchar, prm::Cuchar)

Scales all coefficients by order. If inv == 0, scales coefficients by order (derivation), else scales coefficients by 1/order (integration).

Input

  • t – Source complex TPSA
  • inv – Put order up, divide, scale by inv of value of order
  • prm – Parameters flag. If set to 0x0, the scaling excludes the order of the parameters in the monomials. Else, scaling is with total order of monomial

Output

  • r – Destination complex TPSA
GTPSA.mad_ctpsa_set0!Method
mad_ctpsa_set0!(t::Ptr{CTPSA}, a::ComplexF64, b::ComplexF64)

Sets the 0th order coefficient (scalar part of TPSA) according to coef[0] = a*coef[0] + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • a – Scaling of current 0th order value
  • b – Constant added to current 0th order value
GTPSA.mad_ctpsa_set0_r!Method
mad_ctpsa_set0_r!(t::Ptr{CTPSA}, a_re::Cdouble, a_im::Cdouble, b_re::Cdouble, b_im::Cdouble)

Sets the 0th order coefficient (scalar part of TPSA) according to coef[0] = a*coef[0] + b. Does not modify other values in TPSA. Equivalent to mad_ctpsa_set0 but without complex-by-value arguments.

Input

  • t – TPSA
  • a_re – Real part of a
  • a_im – Imaginary part of a
  • b_re – Real part of b
  • b_im – Imaginary part of b
GTPSA.mad_ctpsa_seti!Method
mad_ctpsa_seti!(t::Ptr{CTPSA}, i::Cint, a::ComplexF64, b::ComplexF64)

Sets the coefficient of monomial at index i to coef[i] = a*coef[i] + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • i – Index of monomial
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_ctpsa_seti_r!Method
mad_ctpsa_seti_r!(t::Ptr{CTPSA}, i::Cint, a_re::Cdouble, a_im::Cdouble, b_re::Cdouble, b_im::Cdouble)

Sets the coefficient of monomial at index i to coef[i] = a*coef[i] + b. Does not modify other values in TPSA. Equivalent to mad_ctpsa_seti but without complex-by-value arguments.

Input

  • t – TPSA
  • i – Index of monomial
  • a_re – Real part of a
  • a_im – Imaginary part of a
  • b_re – Real part of b
  • b_im – Imaginary part of b
GTPSA.mad_ctpsa_setm!Method
mad_ctpsa_setm!(t::Ptr{CTPSA}, n::Cint, m::Vector{Cuchar}, a::ComplexF64, b::ComplexF64)

Sets the coefficient of monomial defined by byte array m to coef = a*coef + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as byte array
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_ctpsa_setm_r!Method
mad_ctpsa_setm_r!(t::Ptr{CTPSA}, n::Cint, m::Vector{Cuchar}, a_re::Cdouble, a_im::Cdouble, b_re::Cdouble, b_im::Cdouble)

Sets the coefficient of monomial defined by byte array m to coef = a*coef + b. Does not modify other values in TPSA. Equivalent to mad_ctpsa_setm but without complex-by-value arguments.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as byte array
  • a_re – Real part of a
  • a_im – Imaginary part of a
  • b_re – Real part of b
  • b_im – Imaginary part of b
GTPSA.mad_ctpsa_setnam!Method
mad_ctpsa_setnam!(t::Ptr{CTPSA}, nam::Cstring)

Sets the name of the CTPSA.

Input

  • t – Complex TPSA
  • nam – Name to set for CTPSA
GTPSA.mad_ctpsa_sets!Method
mad_ctpsa_sets!(t::Ptr{CTPSA}, n::Cint, s::Cstring, a::ComplexF64, b::ComplexF64)

Sets the coefficient of monomial defined by string s to coef = a*coef + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as string
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_ctpsa_sets_r!Method
mad_ctpsa_sets_r!(t::Ptr{CTPSA}, n::Cint, s::Cstring, a_re::Cdouble, a_im::Cdouble, b_re::Cdouble, b_im::Cdouble)

Sets the coefficient of monomial defined by string s to coef = a*coef + b. Does not modify other values in TPSA. Equivalent to mad_ctpsa_set but without complex-by-value arguments.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as string
  • a_re – Real part of a
  • a_im – Imaginary part of a
  • b_re – Real part of b
  • b_im – Imaginary part of b
GTPSA.mad_ctpsa_setsm!Method
mad_ctpsa_setsm!(t::Ptr{CTPSA}, n::Cint, m::Vector{Cint}, a::ComplexF64, b::ComplexF64)

Sets the coefficient of monomial defined by sparse monomial m to coef = a*coef + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as sparse monomial
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_ctpsa_setsm_r!Method
mad_ctpsa_setsm_r!(t::Ptr{CTPSA}, n::Cint, m::Vector{Cint}, a_re::Cdouble, a_im::Cdouble, b_re::Cdouble, b_im::Cdouble)

Sets the coefficient of monomial defined by sparse monomial m to coef = a*coef + b. Does not modify other values in TPSA. Equivalent to mad_ctpsa_setsm but without complex-by-value arguments.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as sparse monomial
  • a_re – Real part of a
  • a_im – Imaginary part of a
  • b_re – Real part of b
  • b_im – Imaginary part of b
GTPSA.mad_ctpsa_setv!Method
mad_ctpsa_setv!(t::Ptr{CTPSA}, i::Cint, n::Cint, v::Vector{ComplexF64})

Vectorized setter of the coefficients for monomials with indices i..i+n. Useful for putting a matrix into a map.

Input

  • t – TPSA
  • i – Starting index of monomials to set coefficients
  • n – Number of monomials to set coefficients of starting at i
  • v – Array of coefficients for monomials i..i+n
GTPSA.mad_ctpsa_setval!Method
mad_ctpsa_setval!(t::Ptr{CTPSA}, v::ComplexF64)

Sets the scalar part of the TPSA to v and all other values to 0 (sets the TPSA order to 0).

Input

  • t – TPSA to set to scalar
  • v – Scalar value to set TPSA
GTPSA.mad_ctpsa_setval_r!Method
mad_ctpsa_setval_r!(t::Ptr{CTPSA}, v_re::Cdouble, v_im::Cdouble)

Sets the scalar part of the TPSA to v and all other values to 0 (sets the TPSA order to 0). Equivalent to mad_ctpsa_setval but without complex-by-value arguments.

Input

  • t – TPSA to set to scalar
  • v_re – Real part of scalar value to set TPSA
  • v_im – Imaginary part of scalar value to set TPSA
GTPSA.mad_ctpsa_setvar!Method
mad_ctpsa_setvar!(t::Ptr{CTPSA}, v::ComplexF64, iv_::Cint, scl_::ComplexF64)

Sets the 0th and 1st order values for the variables.

Input

  • t – Real TPSA
  • v – 0th order value (coefficient)
  • iv_ – (Optional) Variable index, optional if order of TPSA is 0 (behaves like mad_ctpsa_setval then)
  • scl_ – 1st order variable value (typically will be 1)
GTPSA.mad_ctpsa_setvar_r!Method
mad_ctpsa_setvar_r!(t::Ptr{CTPSA}, v_re::Cdouble, v_im::Cdouble, iv_::Cint, scl_re_::Cdouble, scl_im_::Cdouble)

Sets the 0th and 1st order values for the variables. Equivalent to mad_ctpsa_setvar but without complex-by-value arguments.

Input

  • t – Complex TPSA
  • v_re – Real part of 0th order value
  • v_im – Imaginary part of 0th order value
  • iv_ – (Optional) Variable index, optional if order of TPSA is 0 (behaves like mad_ctpsa_setval then)
  • scl_re_ – (Optional) Real part of 1st order variable value
  • scl_im_ – (Optional)Imaginary part of 1st order variable value
GTPSA.mad_ctpsa_sin!Method
mad_ctpsa_sin!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the sin of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sin(a)
GTPSA.mad_ctpsa_sinc!Method
mad_ctpsa_sinc!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the sinc of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sinc(a)
GTPSA.mad_ctpsa_sincos!Method
mad_ctpsa_sincos!(a::Ptr{CTPSA}, s::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA s = sin(a) and TPSA c = cos(a)

Input

  • a – Source TPSA a

Output

  • s – Destination TPSA s = sin(a)
  • c – Destination TPSA c = cos(a)
GTPSA.mad_ctpsa_sincosh!Method
mad_ctpsa_sincosh!(a::Ptr{CTPSA}, s::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA s = sinh(a) and TPSA c = cosh(a)

Input

  • a – Source TPSA a

Output

  • s – Destination TPSA s = sinh(a)
  • c – Destination TPSA c = cosh(a)
GTPSA.mad_ctpsa_sinh!Method
mad_ctpsa_sinh!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the sinh of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sinh(a)
GTPSA.mad_ctpsa_sinhc!Method
mad_ctpsa_sinhc!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the sinhc of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sinhc(a)
GTPSA.mad_ctpsa_sqrt!Method
mad_ctpsa_sqrt!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the sqrt of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sqrt(a)
GTPSA.mad_ctpsa_sub!Method
mad_ctpsa_sub!(a::Ptr{CTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination TPSA c = a - b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a - b
GTPSA.mad_ctpsa_subt!Method
mad_ctpsa_subt!(a::Ptr{CTPSA}, b::Ptr{RTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a - b (internal real-to-complex conversion).

Input

  • a – Source CTPSA a
  • b – Source RTPSA b

Output

  • c – Destination CTPSA c = a - b
GTPSA.mad_ctpsa_tan!Method
mad_ctpsa_tan!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the tan of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = tan(a)
GTPSA.mad_ctpsa_tanh!Method
mad_ctpsa_tanh!(a::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets TPSA c to the tanh of TPSA a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = tanh(a)
GTPSA.mad_ctpsa_taylor!Method
mad_ctpsa_taylor!(a::Ptr{CTPSA}, n::Cint, coef::Vector{ComplexF64}, c::Ptr{CTPSA})

Computes the result of the Taylor series up to order n-1 with Taylor coefficients coef for the scalar value in a. That is, c = coef[0] + coef[1]*a_0 + coef[2]*a_0^2 + ... where a_0 is the scalar part of TPSA a

Input

  • a – TPSA a
  • nOrder-1 of Taylor expansion, size of coef array
  • coef – Array of coefficients in Taylor s
  • c – Result
GTPSA.mad_ctpsa_tdif!Method
mad_ctpsa_tdif!(a::Ptr{RTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

For each homogeneous polynomial in RTPSA a and CTPSA b, calculates either the relative error or absolute error for each order. If the maximum coefficient for a given order in a is > 1, the relative error is computed for that order. Else, the absolute error is computed. This is very useful for comparing maps between codes or doing unit tests. In Julia, essentially:

c_i = (a_i.-b_i)/maximum([abs.(a_i)...,1]) where a_i and b_i are vectors of the monomials for an order i

Input

  • a – Source RTPSA a
  • b – Source CTPSA b

Output

  • c – Destination CTPSA c
GTPSA.mad_ctpsa_tdiv!Method
mad_ctpsa_tdiv!(a::Ptr{RTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a / b (internal real-to-complex conversion).

Input

  • a – Source RTPSA a
  • b – Source CTPSA b

Output

  • c – Destination CTPSA c = a / b
GTPSA.mad_ctpsa_tpoisbra!Method
mad_ctpsa_tpoisbra!(a::Ptr{RTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA}, nv::Cint)

Sets TPSA c to the poisson bracket of RTPSA a and CTPSA b (internal real-to-complex conversion).

Input

  • a – Source RTPSA a
  • b – Source CTPSA b
  • nv – Number of variables in the TPSA

Output

  • c – Destination CTPSA c
GTPSA.mad_ctpsa_tpow!Method
mad_ctpsa_tpow!(a::Ptr{RTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a ^ b (internal real-to-complex conversion).

Input

  • a – Source RTPSA a
  • b – Source CTPSA b

Output

  • c – Destination TPSA c = a ^ b
GTPSA.mad_ctpsa_translate!Method
mad_ctpsa_translate!(na::Cint, ma::Vector{Ptr{CTPSA}}, nb::Cint, tb::Vector{ComplexF64}, mc::Vector{Ptr{CTPSA}})

Translates the expansion point of the map by the amount tb.

Input

  • na – Number of TPSAS in the map
  • ma – Map ma
  • nb – Length of tb
  • tb – Vector of amount to translate for each variable

Output

  • mc – Map evaluated at the new point translated tb from the original evaluation point
GTPSA.mad_ctpsa_tsub!Method
mad_ctpsa_tsub!(a::Ptr{RTPSA}, b::Ptr{CTPSA}, c::Ptr{CTPSA})

Sets the destination CTPSA c = a - b (internal real-to-complex conversion).

Input

  • a – Source RTPSA a
  • b – Source CTPSA b

Output

  • c – Destination CTPSA c = a - b
GTPSA.mad_ctpsa_uid!Method
mad_ctpsa_uid!(t::Ptr{CTPSA}, uid_::Cint)::Cint

Sets the TPSA uid if uid_ != 0, and returns the current (previous if set) TPSA uid.

Input

  • t – Complex TPSA
  • uid_uid to set in the TPSA if uid_ != 0

Output

  • ret – Current (previous if set) TPSA uid
GTPSA.mad_ctpsa_unit!Method
mad_ctpsa_unit!(t::Ptr{CTPSA}, r::Ptr{CTPSA})

Interpreting TPSA a vector, gets the "unit vector", e.g. r = t/norm(t). May be useful for checking for convergence.

Input

  • t – Source TPSA x

Output

  • r – Destination TPSA r
GTPSA.mad_ctpsa_vec2fld!Method
mad_ctpsa_vec2fld!(na::Cint, a::Ptr{CTPSA}, mc::Vector{Ptr{CTPSA}})

Assuming the variables in the TPSA are canonically-conjugate, and ordered so that the canonically- conjugate variables are consecutive (q1, p1, q2, p2, ...), calculates the vector field (Hamilton's equations) from the passed Hamiltonian, defined as [da/dp1, -da/dq1, ...]

Input

  • na – Number of TPSA in mc consistent with number of variables in a
  • a – Hamiltonian as a TPSA

Output

  • mc – Vector field derived from a using Hamilton's equations
GTPSA.mad_desc_del!Method
mad_desc_del!(d_::Ptr{Desc})

Calls the destructor for the passed descriptor.

GTPSA.mad_desc_getnv!Method
mad_desc_getnv!(d::Ptr{Desc}, mo_::Ref{Cuchar}, np_::Ref{Cint}, po_::Ref{Cuchar}::Cint

Returns the number of variables in the descriptor, and sets the passed mo_, np_, and po_ to the maximum order, number of parameters, and parameter order respectively.

Input

  • d – Descriptor

Output

  • mo_ – (Optional) Maximum order of the descriptor
  • np_ – (Optional) Number of parameters of the descriptor
  • po_ – (Optional) Parameter order of the descriptor
  • ret – Number of variables in TPSA
GTPSA.mad_desc_gtrunc!Method
mad_desc_gtrunc!(d::Ptr{Desc}, to::Cuchar)::Cuchar

Sets the global truncation order to of the TPSA, and returns the old global truncation order.

Input

  • d – Descriptor
  • to – New global truncation order

Output

  • oldto – Old global truncation order
GTPSA.mad_desc_idxmMethod
mad_desc_idxm(d::Ptr{Desc}, n::Cint, m::Vector{Cuchar})::Cint

Returns the index of the monomial as byte array m in the descriptor, or -1 if the monomial is invalid.

Input

  • d – Descriptor
  • n – Monomial length
  • m – Monomial as byte array

Output

  • ret – Monomial index or -1 if invalid
GTPSA.mad_desc_idxsMethod
mad_desc_idxs(d::Ptr{Desc}, n::Cint, s::Cstring)::Cint

Returns the index of the monomial as string s in the descriptor, or -1 if the monomial is invalid.

Input

  • d – Descriptor
  • n – String length or 0 if unknown
  • s – Monomial as string "[0-9]*"

Output

  • ret – Monomial index or -1 if invalid monomial
GTPSA.mad_desc_idxsmMethod
mad_desc_idxsm(d::Ptr{Desc}, n::Cint, m::Vector{Cint})::Cint

Returns the index of the monomial as sparse monomial m, indexed as [(i,o)], in the descriptor, or -1 if the monomial is invalid.

Input

  • d – Descriptor
  • n – Monomial length
  • m – Sparse monomial [(idx,ord)]

Output

  • ret – Monomial index or -1 if invalid
GTPSA.mad_desc_infoMethod
mad_desc_info(d::Ptr{Desc}, fp::Ptr{Cvoid})

For debugging.

Input

  • d – Descriptor to debug
  • fp – File to write to. If null, will write to stdout
GTPSA.mad_desc_isvalidmMethod
mad_desc_isvalidm(d::Ptr{Desc}, n::Cint, m::Vector{Cuchar})::Cuchar

Checks if monomial as byte array m is valid given maximum order of descriptor.

Input

  • d – Descriptor
  • n – Length of monomial
  • m – Monomial as byte array

Output

  • ret – True if valid, false if invalid
GTPSA.mad_desc_isvalidsMethod
mad_desc_isvalids(d::Ptr{Desc}, n::Cint, s::Cstring)::Cuchar

Checks if monomial as string s is valid given maximum order of descriptor.

Input

  • d – Descriptor
  • n – Monomial string length
  • s – Monomial as string "[0-9]*"

Output

  • ret – True if valid, false if invalid
GTPSA.mad_desc_isvalidsmMethod
mad_desc_isvalidsm(d::Ptr{Desc}, n::Cint, m::Vector{Cint})::Cuchar

Checks the monomial as sparse monomial m (monomial stored as sequence of integers with each pair [(i,o)] such that i = index, o = order) is valid given the maximum order of the descriptor.

Input

  • d – Descriptor
  • n – Length of monomial
  • m – Sparse monomial [(idx, ord)]

Output

  • ret – True if valid, false if invalid
GTPSA.mad_desc_maxlenMethod
mad_desc_maxlen(d::Ptr{Desc}, mo::Cuchar)::Cint

Gets the maximum length of the TPSA given an order.

Input

  • d – Descriptor
  • mo – Order (ordlen(maxord) == maxlen)

Output

  • ret – monomials in 0..order
GTPSA.mad_desc_maxordMethod
mad_desc_maxord(d::Ptr{Desc}, nn::Cint, no_::Vector{Cuchar})::Cuchar

Sets the order of the variables and parameters of the TPSA to those specified in no_ and returns the maximum order of the TPSA.

Input

  • d – Descriptor
  • nn – Number of variables + number of parameters, no_[1..nn]
  • no_ – (Optional) Orders of parameters to be filled if provided

Output

  • ret – Maximum order of TPSA
GTPSA.mad_desc_mono!Method
mad_desc_mono!(d::Ptr{Desc}, i::Cint, n::Cint, m_::Vector{Cuchar}, p_::Vector{Cuchar})::Cuchar

Returns the order of the monomial at index i, and if n and m_ are provided, then will also fill m_ with the monomial at this index. Also will optionally return the order of the parameters in the monomial if p_ is provided

Input

  • d – Descriptor
  • i – Slot index (must be valid)
  • n – Monomial length (must be provided if m_ is to be filled)

Output

  • ret – Monomial order at slot index
  • m_ – (Optional) Monomial to fill if provided
  • p_ – (Optional) Order of parameters in monomial if provided
GTPSA.mad_desc_newvMethod
mad_desc_newv(nv::Cint, mo::Cuchar)::Ptr{Desc}

Creates a TPSA descriptor with the specified number of variables and maximum order. The number of parameters is set to 0.

Input

  • nv – Number of variables in the TPSA
  • mo – Maximum order of TPSA, mo = max(1, mo)

Output

  • ret – Descriptor with the specified number of variables and maximum order
GTPSA.mad_desc_newvpMethod
mad_desc_newvp(nv::Cint, mo::Cuchar, np_::Cint, po_::Cuchar)::Ptr{Desc}

Creates a TPSA descriptor with the specifed number of variables, maximum order, number of parameters, and parameter order.

Input

  • nv – Number of variables
  • mo – Maximum order of TPSA INCLUDING PARAMETERS, mo = max(1, mo)
  • np_ – (Optional) Number of parameters, default is 0
  • po_ – (Optional) Order of parameters, po = max(1, po_)

Output

  • ret – Descriptor with the specified nv, mo, np, and po
GTPSA.mad_desc_newvpoMethod
mad_desc_newvpo(nv::Cint, mo::Cuchar, np_::Cint, po_::Cuchar, no_::Vector{Cuchar})::Ptr{Desc}

Creates a TPSA descriptor with the specifed number of variables, maximum order for both variables and parameters, number of parameters, parameter order, and individual variable/parameter orders specified in no. The first nv entries in no correspond to the variables' orders and the next np entries correspond the parameters' orders.

Input

  • nv – Number of variables
  • mo – Maximum order of TPSA (mo = max(mo , no[0 :nn-1]), nn = nv+np)
  • np_ – (Optional) Number of parameters, default is 0
  • po_ – (Optional) Order of parameters (po = max(po_, no[nv:nn-1]), po <= mo)
  • no_ – (Optional) Array of orders of variables and parameters

Output

  • ret – Descriptor with the specified nv, mo, np, po, no.
GTPSA.mad_desc_nxtbyordMethod
mad_desc_nxtbyord(d::Ptr{Desc}, n::Cint, m::Vector{Cuchar})::Cint

Returns the next monomial after monomial m in the TPSA when sorted by order.

Input

  • d – Descriptor
  • n – Monomial length
  • m – Monomial as byte array

Output

  • idx – Monomial index or -1 if no valid next monomial
GTPSA.mad_desc_nxtbyvarMethod
mad_desc_nxtbyvar(d::Ptr{Desc}, n::Cint, m::Vector{Cuchar})::Cint

Returns the next monomial after monomial m in the TPSA when sorted by variable.

Input

  • d – Descriptor
  • n – Monomial length
  • m – Monomial as byte array

Output

  • idx – Monomial index or -1 if no valid next monomial
GTPSA.mad_mono_add!Method
mad_mono_add!(n::Cint, a::Vector{Cuchar}, b::Vector{Cuchar}, r::Vector{Cuchar})

Sets monomial r = a + b.

Input

  • n – Length of monomials
  • a – Source monomial a
  • b – Source monomial b

Output

  • r – Destination monomial, r = a + b
GTPSA.mad_mono_cat!Method
mad_mono_cat!(n::Cint, a::Vector{Cuchar}, m::Cint, b::Vector{Cuchar}, r::Vector{Cuchar})

Sets monomial r equal to the concatenation of the monomials a and b

Input

  • n – Length of monomonial a
  • a – Source monomial a
  • m – Length of monomial b
  • b – Source monomial b

Output

  • r – Destination monomial of concatenation of a and b (length n+m)
GTPSA.mad_mono_cmpMethod
mad_mono_cmp(n::Cint, a::Vector{Cuchar}, b::Vector{Cuchar})::Cint

Compares monomial a to monomial b, and returns the first difference in the lowest order variables.

Input

  • n – Length of monomials
  • a – Monomial a
  • b – Monomial b

Output

  • ret – First a[i]-b[i] != 0
GTPSA.mad_mono_copy!Method
mad_mono_copy!(n::Cint, a::Vector{Cuchar}, r::Vector{Cuchar})

Copies monomial a to monomial r.

Input

  • n – Length of monomials
  • a – Source monomial
  • r – Destination monomial
GTPSA.mad_mono_eqMethod
mad_mono_eq(n::Cint, a::Vector{Cuchar}, b::Vector{Cuchar})::Cuchar

Checks if the monomial a is equal to the monomial b.

Input

  • n – Length of monomials
  • a – Monomial a
  • b – Monomial b

Output

  • ret – True if the monomials are equal, false if otherwise
GTPSA.mad_mono_fill!Method
mad_mono_fill!(n::Cint, a::Vector{Cuchar}, v::Cuchar)

Fills the monomial a with the value v.

Input

  • n – Monomial length
  • a – Monomial
  • v – Value
GTPSA.mad_mono_leMethod
mad_mono_le(n::Cint, a::Vector{Cuchar}, b::Vector{Cuchar})::Cuchar

Checks if monomial a is less than or equal to monomial b.

Input

  • n – Length of monomials
  • a – Monomial a
  • b – Monomial b

Output

  • ret – True if a <= mono_b, false otherwise
GTPSA.mad_mono_ltMethod
mad_mono_lt(n::Cint, a::Vector{Cuchar}, b::Vector{Cuchar})::Cuchar

Checks if monomial a is less than monomial b.

Input

  • n – Length of monomials
  • a – Monomial a
  • b – Monomial b

Output

  • ret – True if a < b, false otherwise
GTPSA.mad_mono_maxMethod
mad_mono_max(n::Cint, a::Vector{Cuchar})::Cuchar

Returns the maximum order of the monomial.

Input

  • n – Length of monomial
  • a – Monomial

Output

  • mo – Maximum order of monomial a
GTPSA.mad_mono_minMethod
mad_mono_min(n::Cint, a::Vector{Cuchar})::Cuchar

Returns the minimum order of the monomial.

Input

  • n – Length of monomial
  • a – Monomial

Output

  • mo – Mininum order of monomial a
GTPSA.mad_mono_ordMethod
mad_mono_ord(n::Cint, a::Vector{Cuchar})::Cint

Returns the sum of the orders of the monomial a.

Input

  • n – Monomial length
  • a – Monomial

Output

  • s – Sum of orders of monomial
GTPSA.mad_mono_ordpMethod
mad_mono_ordp(n::Cint, a::Vector{Cuchar}, stp::Cint)::Cdouble

Returns the product of each stp-th order in monomial a. For example, stp = 2 collects every order in the monomial with a step of 2 between each. As a is a pointer, the product can be started at any element in the monomial.

Input

  • n – Monomial length
  • a – Monomial as byte array
  • stp – Step over which orders to include in the product

Output

  • p – Product of orders of monomial separated by stp.
GTPSA.mad_mono_ordpfMethod
mad_mono_ordpf(n::Cint, a::Vector{Cuchar}, stp::Cint)::Cdouble

Returns the product of factorials each stp-th order in monomial a. For example, stp = 2 collects every order in the monomial with a step of 2 between each. As a is a pointer, the product can be started at any element in the monomial.

Input

  • n – Monomial length
  • a – Monomial as byte array
  • stp – Step over which orders to include in the product of factorials

Output

  • p – Product of factorials of orders of monomial separated by stp
GTPSA.mad_mono_printMethod
mad_mono_print(n::Cint, a::Vector{Cuchar}, fp_::Ptr{Cvoid})

Prints the monomial to stdout.

Input

  • n – Length of monomial
  • a – Source monomial to print to stdout
  • fp_ – C FILE pointer, if null will print to stdout
GTPSA.mad_mono_prt!Method
mad_mono_prt(n::Cint, a::Vector{Cuchar}, s::Ptr{Cuchar})::Cstring

Writes the monomial defined by the byte array a (with orders stored as hexadecimal) into a null terminated string s.

Input

  • n – Monomial and string length
  • a – Monomial as byte array

Output

  • ret – Monomial as string
GTPSA.mad_mono_rcmpMethod
mad_mono_rcmp(n::Cint, a::Vector{Cuchar}, b::Vector{Cuchar})::Cint

Compares monomial a to monomial b starting from the right (when the monomials are ordered by variable, which is almost never the case) and returns the first difference in the lowest order variables.

Input

  • n – Length of monomials
  • a – Monomial a
  • b – Monomial b

Output

  • ret – First a[i]-b[i] != 0 where i starts from the end.
GTPSA.mad_mono_rev!Method
mad_mono_rev!(n::Cint, a::Vector{Cuchar}, r::Vector{Cuchar})

Sets destination monomial r equal to the reverse of source monomial a.

Input

  • n – Lengths of monomials
  • a – Source monomial a

Output

  • r – Destination monomial of reverse monomial a
GTPSA.mad_mono_str!Method
mad_mono_str!(n::Cint, a::Vector{Cuchar}, s::Cstring)::Cint

Writes the monomial defined in the string s, which stores the orders in a human-readable format (e.g. 10 is 10, not 0xa), into the byte array a with the orders specified in hexadecimal.

Input

  • n – Monomial and string length
  • s – Monomial as string "[0-9]*"

Output

  • a – Monomial as a byte array converted from the input string
  • i – Adjusted size n of byte array if '' found
GTPSA.mad_mono_sub!Method
mad_mono_sub!(n::Cint, a::Vector{Cuchar}, b::Vector{Cuchar}, r::Vector{Cuchar})

Sets monomial r = a - b.

Input

  • n – Length of monomials
  • a – Source monomial a
  • b – Source monomial b

Output

  • r – Destination monomial, r = a - b
GTPSA.mad_tpsa_abs!Method
mad_tpsa_abs!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the absolute value of TPSA a. Specifically, the result contains a TPSA with the abs of all coefficients.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = |a|
GTPSA.mad_tpsa_acc!Method
mad_tpsa_acc!(a::Ptr{RTPSA}, v::Cdouble, c::Ptr{RTPSA})

Adds a*v to TPSA c. Aliasing OK.

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c += v*a
GTPSA.mad_tpsa_acos!Method
mad_tpsa_acos!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the acos of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = acos(a)
GTPSA.mad_tpsa_acosh!Method
mad_tpsa_acosh!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the acosh of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA `c = acosh(a)'
GTPSA.mad_tpsa_acot!Method
mad_tpsa_acot!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the acot of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = acot(a)
GTPSA.mad_tpsa_acoth!Method
mad_tpsa_acoth!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the acoth of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA `c = acoth(a)'
GTPSA.mad_tpsa_add!Method
mad_tpsa_add!(a::Ptr{RTPSA}, b::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets the destination TPSA c = a + b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a + b
GTPSA.mad_tpsa_asin!Method
mad_tpsa_asin!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the asin of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = asin(a)
GTPSA.mad_tpsa_asinc!Method
mad_tpsa_asinc!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the asinc(a) = asin(a)/a

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = asinc(a) = asin(a)/a
GTPSA.mad_tpsa_asinh!Method
mad_tpsa_asinh!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the asinh of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA `c = asinh(a)'
GTPSA.mad_tpsa_asinhc!Method
mad_tpsa_asinhc!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the asinhc of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA `c = asinhc(a)'
GTPSA.mad_tpsa_atan!Method
mad_tpsa_atan!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the atan of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = atan(a)
GTPSA.mad_tpsa_atan2!Method
mad_tpsa_atan2!(y::Ptr{RTPSA}, x::Ptr{RTPSA}, r::Ptr{RTPSA})

Sets TPSA r to atan2(y,x)

Input

  • y – Source TPSA y
  • x – Source TPSA x

Output

  • r – Destination TPSA r = atan2(y,x)
GTPSA.mad_tpsa_atanh!Method
mad_tpsa_atanh!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the atanh of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA `c = atanh(a)'
GTPSA.mad_tpsa_ax2pby2pcz2!Method
mad_tpsa_ax2pby2pcz2!(a::Cdouble, x::Ptr{RTPSA}, b::Cdouble, y::Ptr{RTPSA}, c::Cdouble, z::Ptr{RTPSA}, r::Ptr{RTPSA})

r = a*x^2 + b*y^2 + c*z^2

Input

  • a – Scalar a
  • x – TPSA x
  • b – Scalar b
  • y – TPSA y
  • c – Scalar c
  • z – TPSA z

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_axpb!Method
mad_tpsa_axpb!(a::Cdouble, x::Ptr{RTPSA}, b::Cdouble, r::Ptr{RTPSA})

r = a*x + b

Input

  • a – Scalar a
  • x – TPSA x
  • b – Scalar b

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_axpbypc!Method
mad_tpsa_axpbypc!(a::Cdouble, x::Ptr{RTPSA}, b::Cdouble, y::Ptr{RTPSA}, c::Cdouble, r::Ptr{RTPSA})

r = a*x + b*y + c

Input

  • a – Scalar a
  • x – TPSA x
  • b – Scalar b
  • y – TPSA y
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_axpsqrtbpcx2!Method
mad_tpsa_axpsqrtbpcx2!(x::Ptr{RTPSA}, a::Cdouble, b::Cdouble, c::Cdouble, r::Ptr{RTPSA})

r = a*x + sqrt(b + c*x^2)

Input

  • x – TPSA x
  • a – Scalar a
  • b – Scalar b
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_axypb!Method
mad_tpsa_axypb!(a::Cdouble, x::Ptr{RTPSA}, y::Ptr{RTPSA}, b::Cdouble, r::Ptr{RTPSA})

r = a*x*y + b

Input

  • a – Scalar a
  • x – TPSA x
  • y – TPSA y
  • b – Scalar b

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_axypbvwpc!Method
mad_tpsa_axypbvwpc!(a::Cdouble, x::Ptr{RTPSA}, y::Ptr{RTPSA}, b::Cdouble, v::Ptr{RTPSA}, w::Ptr{RTPSA}, c::Cdouble, r::Ptr{RTPSA})

r = a*x*y + b*v*w + c

Input

  • a – Scalar a
  • x – TPSA x
  • y – TPSA y
  • b – Scalar b
  • v – TPSA v
  • w – TPSA w
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_axypbzpc!Method
mad_tpsa_axypbzpc!(a::Cdouble, x::Ptr{RTPSA}, y::Ptr{RTPSA}, b::Cdouble, z::Ptr{RTPSA}, c::Cdouble, r::Ptr{RTPSA})

r = a*x*y + b*z + c

Input

  • a – Scalar a
  • x – TPSA x
  • y – TPSA y
  • b – Scalar b
  • z – TPSA z
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_clear!Method
mad_tpsa_clear!(t::Ptr{RTPSA})

Clears the TPSA (reset to 0)

Input

  • t – TPSA
GTPSA.mad_tpsa_compose!Method
mad_tpsa_compose!(na::Cint, ma::Vector{Ptr{RTPSA}}, nb::Cint, mb::Vector{Ptr{RTPSA}}, mc::Vector{Ptr{RTPSA}})

Composes two maps.

Input

  • na – Number of TPSAs in map ma
  • ma – map ma
  • nb – Number of TPSAs in map mb
  • mb – map mb

Output

  • mc – Composition of maps ma and mb
GTPSA.mad_tpsa_convert!Method
mad_tpsa_convert!(t::Ptr{RTPSA}, r::Ptr{RTPSA}, n::Cint, t2r_::Vector{Cint}, pb::Cint)

General function to convert TPSAs to different orders and reshuffle canonical coordinates. The destination TPSA will be of order n, and optionally have the variable reshuffling defined by t2r_ and poisson bracket sign. e.g. if t2r_ = {1,2,3,4,6,5} and pb = -1, canonical coordinates 6 and 5 are swapped and the new 5th canonical coordinate will be negated. Useful for comparing with different differential algebra packages.

Input

  • t – Source TPSA
  • n – Length of vector
  • t2r_ – (Optional) Vector of index lookup
  • pb – Poisson bracket, 0, 1:fwd, -1:bwd

Output

  • r – Destination TPSA with specified order and canonical coordinate reshuffling.
GTPSA.mad_tpsa_copy!Method
mad_tpsa_copy!(t::Ptr{RTPSA}, r::Ptr{RTPSA})

Makes a copy of the TPSA t to r.

Input

  • t – Source TPSA

Output

  • r – Destination TPSA
GTPSA.mad_tpsa_cos!Method
mad_tpsa_cos!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the cos of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = cos(a)
GTPSA.mad_tpsa_cosh!Method
mad_tpsa_cosh!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the cosh of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = cosh(a)
GTPSA.mad_tpsa_cot!Method
mad_tpsa_cot!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the cot of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = cot(a)
GTPSA.mad_tpsa_coth!Method
mad_tpsa_coth!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the coth of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = coth(a)
GTPSA.mad_tpsa_cutord!Method
mad_tpsa_cutord!(t::Ptr{RTPSA}, r::Ptr{RTPSA}, ord::Cint)

Cuts the TPSA off at the given order and above, or if ord is negative, will cut orders below abs(ord) (e.g. if ord = -3, then orders 0-3 are cut off).

Input

  • t – Source TPSA
  • ord – Cut order: 0..-ord or ord..mo

Output

  • r – Destination TPSA
GTPSA.mad_tpsa_cycle!Method
mad_tpsa_cycle!(t::Ptr{RTPSA}, i::Cint, n::Cint, m_::Vector{Cuchar}, v_::Ref{Cdouble})::Cint

Used for scanning through each nonzero monomial in the TPSA. Given a starting index (-1 if starting at 0), will optionally fill monomial m_ with the monomial at index i and the value at v_ with the monomials coefficient, and return the next NONZERO monomial index in the TPSA. This is useful for building an iterator through the TPSA.

Input

  • t – TPSA to scan
  • i – Index to start from (-1 to start at 0)
  • n – Length of monomial
  • m_ – (Optional) Monomial to be filled if provided
  • v_ – (Optional) Pointer to value of coefficient

Output

  • i – Index of next nonzero monomial in the TPSA, or -1 if reached the end
GTPSA.mad_tpsa_debugMethod
mad_tpsa_debug(t::Ptr{RTPSA}, name_::Cstring, fnam_::Cstring, line_::Cint, stream_::Ptr{Cvoid})

Prints TPSA with all information of data structure.

Input

  • t – TPSA
  • name_ – (Optional) Name of TPSA
  • fnam_ – (Optional) File name to print to
  • line_ – (Optional) Line number in file to start at
  • stream_ – (Optional) I/O stream to print to, default is stdout
GTPSA.mad_tpsa_del!Method
mad_tpsa_del!(t::Ptr{RTPSA})

Calls the destructor for the TPSA.

Input

  • t – TPSA to destruct
GTPSA.mad_tpsa_deriv!Method
mad_tpsa_deriv!(a::Ptr{RTPSA}, c::Ptr{RTPSA}, iv::Cint)

Differentiates TPSA with respect to the variable with index iv.

Input

  • a – Source TPSA to differentiate
  • iv – Index of variable to take derivative wrt to (e.g. derivative wrt x, iv = 1).

Output

  • c – Destination TPSA
GTPSA.mad_tpsa_derivm!Method
mad_tpsa_derivm!(a::Ptr{RTPSA}, c::Ptr{RTPSA}, n::Cint, m::Vector{Cuchar})

Differentiates TPSA with respect to the monomial defined by byte array m.

Input

  • a – Source TPSA to differentiate
  • n – Length of monomial to differentiate wrt
  • m – Monomial to take derivative wrt

Output

  • c – Destination TPSA
GTPSA.mad_tpsa_descMethod
mad_tpsa_desc(t::Ptr{RTPSA})::Ptr{Desc}

Gets the descriptor for the TPSA.

Input

  • t – TPSA

Output

  • ret – Descriptor for the RTPSA
GTPSA.mad_tpsa_dif!Method
mad_tpsa_dif!(a::Ptr{RTPSA}, b::Ptr{RTPSA}, c::Ptr{RTPSA})

For each homogeneous polynomial in TPSAs a and b, calculates either the relative error or absolute error for each order. If the maximum coefficient for a given order in a is > 1, the relative error is computed for that order. Else, the absolute error is computed. This is very useful for comparing maps between codes or doing unit tests. In Julia, essentially:

c_i = (a_i.-b_i)/maximum([abs.(a_i)...,1]) where a_i and b_i are vectors of the monomials for an order i

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c
GTPSA.mad_tpsa_div!Method
mad_tpsa_div!(a::Ptr{RTPSA}, b::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets the destination TPSA c = a / b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a / b
GTPSA.mad_tpsa_equMethod
mad_tpsa_equ(a::Ptr{RTPSA}, b::Ptr{RTPSA}, tol_::Cdouble)::Cuchar

Checks if the TPSAs a and b are equal within the specified tolerance tol_. If tol_ is not specified, DBL_GTPSA.show_epsILON is used.

Input

  • a – TPSA a
  • b – TPSA b
  • tol_ – (Optional) Difference below which the TPSAs are considered equal

Output

  • ret - True if a == b within tol_
GTPSA.mad_tpsa_erf!Method
mad_tpsa_erf!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the erf of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA `c = erf(a)'
GTPSA.mad_tpsa_erfc!Method
mad_tpsa_erfc!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the erfc of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA `c = erfc(a)'
GTPSA.mad_tpsa_eval!Method
mad_tpsa_eval!(na::Cint, ma::Vector{Ptr{RTPSA}}, nb::Cint, tb::Vector{Cdouble}, tc::Vector{Cdouble})

Evaluates the map at the point tb

Input

  • na – Number of TPSAs in the map
  • ma – map ma
  • nb – Length of tb
  • tb – Point at which to evaluate the map

Output

  • tc – Values for each TPSA in the map evaluated at the point tb
GTPSA.mad_tpsa_exp!Method
mad_tpsa_exp!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the exponential of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = exp(a)
GTPSA.mad_tpsa_exppb!Method
mad_tpsa_exppb!(na::Cint, ma::Vector{Ptr{RTPSA}}, mb::Vector{Ptr{RTPSA}}, mc::Vector{Ptr{RTPSA}})

Computes the exponential of fgrad of the vector fields ma and mb, literally exppb(ma, mb) = mb + fgrad(ma, mb) + fgrad(ma, fgrad(ma, mb))/2! + ...

Input

  • na – Length of ma and mb
  • ma – Vector of TPSA ma
  • mb – Vector of TPSA mb

Output

  • mc – Destination vector of TPSA mc
GTPSA.mad_tpsa_fgrad!Method
mad_tpsa_fgrad!(na::Cint, ma::Vector{Ptr{RTPSA}}, b::Ptr{RTPSA}, c::Ptr{RTPSA})

Calculates dot(ma, grad(b))

Input

  • na – Length of ma consistent with number of variables in b
  • ma – Vector of TPSA
  • b – TPSA

Output

  • cdot(ma, grad(b))
GTPSA.mad_tpsa_fld2vec!Method
mad_tpsa_fld2vec!(na::Cint, ma::Vector{Ptr{RTPSA}}, c::Ptr{RTPSA})

Assuming the variables in the TPSA are canonically-conjugate, and ordered so that the canonically- conjugate variables are consecutive (q1, p1, q2, p2, ...), calculates the Hamiltonian one obtains from ther vector field (in the form [da/dp1, -da/dq1, ...])

Input

  • na – Number of TPSA in ma consistent with number of variables in c
  • ma – Vector field

Output

  • c – Hamiltonian as a TPSA derived from the vector field ma
GTPSA.mad_tpsa_get0Method
mad_tpsa_get0(t::Ptr{RTPSA})::Cdouble

Gets the 0th order (scalar) value of the TPSA

Input

  • t – TPSA

Output

  • ret – Scalar value of TPSA
GTPSA.mad_tpsa_getiMethod
mad_tpsa_geti(t::Ptr{RTPSA}, i::Cint)::Cdouble

Gets the coefficient of the monomial at index i. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • i – Monomial index

Output

  • ret – Coefficient of monomial at index i
GTPSA.mad_tpsa_getmMethod
mad_tpsa_getm(t::Ptr{RTPSA}, n::Cint, m::Vector{Cuchar})::Cdouble

Gets the coefficient of the monomial m defined as a byte array. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as byte array

Output

  • ret – Coefficient of monomial m in TPSA
GTPSA.mad_tpsa_getord!Method
mad_tpsa_getord!(t::Ptr{RTPSA}, r::Ptr{RTPSA}, ord::Cuchar)

Extract one homogeneous polynomial of the given order

Input

  • t – Source TPSA
  • ord – Order to retrieve

Output

  • r – Destination TPSA
GTPSA.mad_tpsa_getsMethod
mad_tpsa_gets(t::Ptr{RTPSA}, n::Cint, s::Cstring)::Cdouble

Gets the coefficient of the monomial s defined as a string. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as string

Output

  • ret – Coefficient of monomial s in TPSA
GTPSA.mad_tpsa_getsmMethod
mad_tpsa_getsm(t::Ptr{RTPSA}, n::Cint, m::Vector{Cint})::Cdouble

Gets the coefficient of the monomial m defined as a sparse monomial. Generally should use mad_tpsa_cycle instead of this.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as sparse monomial

Output

  • ret – Coefficient of monomial m in TPSA
GTPSA.mad_tpsa_getv!Method
mad_tpsa_getv!(t::Ptr{RTPSA}, i::Cint, n::Cint, v::Vector{Cdouble})

Vectorized getter of the coefficients for monomials with indices i..i+n. Useful for extracting the 1st order parts of a TPSA to construct a matrix (i = 1, n = nv+np = nn).

Input

  • t – TPSA
  • i – Starting index of monomials to get coefficients
  • n – Number of monomials to get coefficients of starting at i

Output

  • v – Array of coefficients for monomials i..i+n
GTPSA.mad_tpsa_hypot!Method
mad_tpsa_hypot!(x::Ptr{RTPSA}, y::Ptr{RTPSA}, r::Ptr{RTPSA})

Sets TPSA r to sqrt(x^2+y^2). Used to oversimplify polymorphism in code but not optimized

Input

  • x – Source TPSA x
  • y – Source TPSA y

Output

  • r – Destination TPSA r = sqrt(x^2+y^2)
GTPSA.mad_tpsa_hypot3!Method
mad_tpsa_hypot3!(x::Ptr{RTPSA}, y::Ptr{RTPSA}, z::Ptr{RTPSA}, r::Ptr{RTPSA})

Sets TPSA r to sqrt(x^2+y^2+z^2). Does NOT allow for r = x, y, z !!!

Input

  • x – Source TPSA x
  • y – Source TPSA y
  • z – Source TPSA z

Output

  • r – Destination TPSA r = sqrt(x^2+y^2+z^2)
GTPSA.mad_tpsa_idxmMethod
mad_tpsa_idxm(t::Ptr{RTPSA}, n::Cint, m::Vector{Cuchar})::Cint

Returns index of monomial in the TPSA given the monomial as a byte array

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as byte array

Output

  • ret – Index of monomial in TPSA
GTPSA.mad_tpsa_idxsMethod
mad_tpsa_idxs(t::Ptr{RTPSA}, n::Cint, s::Cstring)::Cint

Returns index of monomial in the TPSA given the monomial as string. This generally should not be used, as there are no assumptions about which monomial is attached to which index.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as string

Output

  • ret – Index of monomial in TPSA
GTPSA.mad_tpsa_idxsmMethod
mad_tpsa_idxsm(t::Ptr{RTPSA}, n::Cint, m::Vector{Cint})::Cint

Returns index of monomial in the TPSA given the monomial as a sparse monomial. This generally should not be used, as there are no assumptions about which monomial is attached to which index.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as sparse monomial

Output

  • ret – Index of monomial in TPSA
GTPSA.mad_tpsa_init!Method
mad_tpsa_init(t::Ptr{RTPSA}, d::Ptr{Desc}, mo::Cuchar)::Ptr{RTPSA}

Unsafe initialization of an already existing TPSA t with maximum order mo to the descriptor d. mo must be less than the maximum order of the descriptor. t is modified in place and also returned.

Input

  • t – TPSA to initialize to descriptor d
  • d – Descriptor
  • mo – Maximum order of the TPSA (must be less than maximum order of the descriptor)

Output

  • t – TPSA initialized to descriptor d with maximum order mo
GTPSA.mad_tpsa_integ!Method
mad_tpsa_integ!(a::Ptr{RTPSA}, c::Ptr{RTPSA}, iv::Cint)

Integrates TPSA with respect to the variable with index iv.

Input

  • a – Source TPSA to integrate
  • iv – Index of variable to integrate over (e.g. integrate over x, iv = 1).

Output

  • c – Destination TPSA
GTPSA.mad_tpsa_inv!Method
mad_tpsa_inv!(a::Ptr{RTPSA},  v::Cdouble, c::Ptr{RTPSA})

Sets TPSA c to v/a.

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c = v/a
GTPSA.mad_tpsa_invsqrt!Method
mad_tpsa_invsqrt!(a::Ptr{RTPSA}, v::Cdouble, c::Ptr{RTPSA})

Sets TPSA c to v/sqrt(a).

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c = v/sqrt(a)
GTPSA.mad_tpsa_isnulMethod
mad_tpsa_isnul(t::Ptr{RTPSA})::Cuchar

Checks if TPSA is 0 or not

Input

  • t – TPSA to check

Output

  • ret – True or false
GTPSA.mad_tpsa_isvalidMethod
mad_tpsa_isvalid(t::Ptr{RTPSA})::Cuchar

Sanity check of the TPSA integrity.

Input

  • t – TPSA to check if valid

Output

  • ret – True if valid TPSA, false otherwise
GTPSA.mad_tpsa_lenMethod
mad_tpsa_len(t::Ptr{RTPSA})::Cint

Gets the length of the TPSA itself (e.g. the descriptor may be order 10 but TPSA may only be order 2)

Input

  • t – TPSA

Output

  • ret – Length of RTPSA
GTPSA.mad_tpsa_liebra!Method
mad_tpsa_liebra!(na::Cint, ma::Vector{Ptr{RTPSA}}, mb::Vector{Ptr{RTPSA}}, mc::Vector{Ptr{RTPSA}})

Computes the Lie bracket of the vector fields ma and mb, defined as sumi mai (dmb/dxi) - mbi (dma/dx_i).

Input

  • na – Length of ma and mb
  • ma – Vector of TPSA ma
  • mb – Vector of TPSA mb

Output

  • mc – Destination vector of TPSA mc
GTPSA.mad_tpsa_log!Method
mad_tpsa_log!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the log of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = log(a)
GTPSA.mad_tpsa_logaxpsqrtbpcx2!Method
mad_tpsa_logaxpsqrtbpcx2!(x::Ptr{RTPSA}, a::Cdouble, b::Cdouble, c::Cdouble, r::Ptr{RTPSA})

r = log(a*x + sqrt(b + c*x^2))

Input

  • x – TPSA x
  • a – Scalar a
  • b – Scalar b
  • c – Scalar c

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_logpb!Method
mad_tpsa_logpb!(na::Cint, ma::Vector{Ptr{RTPSA}}, mb::Vector{Ptr{RTPSA}}, mc::Vector{Ptr{RTPSA}})

Computes the log of the Poisson bracket of the vector of TPSA ma and mb; the result is the vector field F used to evolve to ma from mb.

Input

  • na – Length of ma and mb
  • ma – Vector of TPSA ma
  • mb – Vector of TPSA mb

Output

  • mc – Destination vector of TPSA mc
GTPSA.mad_tpsa_logxdy!Method
mad_tpsa_logxdy!(x::Ptr{RTPSA}, y::Ptr{RTPSA}, r::Ptr{RTPSA})

r = log(x / y)

Input

  • x – TPSA x
  • y – TPSA y

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_maxord!Method
mad_tpsa_maxord!(t::Ptr{RTPSA}, n::Cint, idx_::Vector{Cint})::Cint

Returns the index to the monomial with maximum abs(coefficient) in the TPSA for all orders 0 to n. If idx_ is provided, it is filled with the indices for the maximum abs(coefficient) monomial for each order up to n.

Input

  • t – TPSA
  • n – Highest order to include in finding the maximum abs(coefficient) in the TPSA, length of idx_ if provided

Output

  • idx_ – (Optional) If provided, is filled with indices to the monomial for each order up to n with maximum abs(coefficient)
  • mi – Index to the monomial in the TPSA with maximum abs(coefficient)
GTPSA.mad_tpsa_mconv!Method
mad_tpsa_mconv!(na::Cint, ma::Vector{Ptr{RTPSA}}, nc::Cint, mc::Vector{Ptr{RTPSA}}, n::Cint, t2r_::Vector{Cint}, pb::Cint)

Equivalent to mad_tpsa_convert, but applies the conversion to all TPSAs in the map ma.

Input

  • na – Number of TPSAs in the map
  • ma – map ma
  • nc – Number of TPSAs in the output map mc
  • n – Length of vector (size of t2r_)
  • t2r_ – (Optional) Vector of index lookup
  • pb – Poisson bracket, 0, 1:fwd, -1:bwd

Output

  • mc – map mc with specified conversions
GTPSA.mad_tpsa_minv!Method
mad_tpsa_minv!(na::Cint, ma::Vector{Ptr{RTPSA}}, mc::Vector{Ptr{RTPSA}})

Inverts the map.

Input

  • na – Number of TPSAs in the map
  • ma – map ma

Output

  • mc – Inversion of map ma
GTPSA.mad_tpsa_mnrmMethod
mad_tpsa_mnrm(na::Cint, ma::Vector{Ptr{RTPSA}})::Cdouble

Computes the norm of the map (sum of absolute value of coefficients of all TPSAs in the map).

Input

  • na – Number of TPSAs in the map
  • ma – map ma

Output

  • nrm – Norm of map (sum of absolute value of coefficients of all TPSAs in the map)
GTPSA.mad_tpsa_mono!Method
mad_tpsa_mono!(t::Ptr{RTPSA}, i::Cint, n::Cint, m_::Vector{Cuchar}, p_::Vector{Cuchar})::Cuchar

Returns the order of the monomial at index i in the TPSA and optionally the monomial at that index is returned in m_ and the order of parameters in the monomial in p_

Input

  • t – TPSA
  • i – Index valid in TPSA
  • n – Length of monomial

Output

  • m_ – (Optional) Monomial at index i in TPSA
  • p_ – (Optional) Order of parameters in monomial
  • ret – Order of monomial in TPSA at index i
GTPSA.mad_tpsa_mul!Method
mad_tpsa_mul!(a::Ptr{RTPSA}, b::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets the destination TPSA c = a * b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a * b
GTPSA.mad_tpsa_namMethod
mad_tpsa_nam(t::Ptr{RTPSA})::Cstring

Get the name of the TPSA.

Input

  • t – TPSA

Output

  • ret – Name of RTPSA (null terminated in C)
GTPSA.mad_tpsa_newMethod
mad_tpsa_new(t::Ptr{RTPSA}, mo::Cuchar)::Ptr{RTPSA}

Creates a blank TPSA with same number of variables/parameters of the inputted TPSA, with maximum order specified by mo. If MAD_TPSA_SAME is passed for mo, the mo currently in t is used for the created TPSA. Ok with t=(tpsa_t*)ctpsa

Input

  • t – TPSA
  • mo – Maximum order of new TPSA

Output

  • ret – New blank TPSA with maximum order mo
GTPSA.mad_tpsa_newdMethod
mad_tpsa_newd(d::Ptr{Desc}, mo::Cuchar)::Ptr{RTPSA}

Creates a TPSA defined by the specified descriptor and maximum order. If MAD_TPSA_DEFAULT is passed for mo, the mo defined in the descriptor is used. If mo > d_mo, then mo = d_mo.

Input

  • d – Descriptor
  • mo – Maximum order

Output

  • t – New TPSA defined by the descriptor
GTPSA.mad_tpsa_nrmMethod
mad_tpsa_nrm(a::Ptr{RTPSA})::Cdouble

Calculates the 1-norm of TPSA a (sum of abs of all coefficients)

Input

  • a – TPSA

Output

  • nrm – 1-Norm of TPSA
GTPSA.mad_tpsa_ordMethod
mad_tpsa_ord(t::Ptr{RTPSA})::Cuchar

Gets the TPSA order.

Input

  • t – TPSA

Output

  • ret – Order of TPSA
GTPSA.mad_tpsa_ordnMethod
mad_tpsa_ordn(n::Cint, t::Vector{Ptr{RTPSA}})::Cuchar

Returns the max order of all TPSAs in t.

Input

  • n – Number of TPSAs
  • t – Array of TPSAs

Output

  • mo – Maximum order of all TPSAs
GTPSA.mad_tpsa_ordvMethod
mad_tpsa_ordv(t::Ptr{RTPSA}, ts::Ptr{RTPSA}...)::Cuchar

Returns maximum order of all TPSAs provided.

Input

  • t – TPSA
  • ts – Variable number of TPSAs passed as parameters

Output

  • mo – Maximum order of all TPSAs provided
GTPSA.mad_tpsa_pminv!Method
mad_tpsa_pminv!(na::Cint, ma::Vector{Ptr{RTPSA}}, mc::Vector{Ptr{RTPSA}}, select::Vector{Cint})

Computes the partial inverse of the map with only the selected variables, specified by 0s or 1s in select.

Input

  • na – Number of TPSAs in ma
  • ma – map ma
  • select – Array of 0s or 1s defining which variables to do inverse on (atleast same size as na)'

Output

  • mc – Partially inverted map using variables specified as 1 in the select array
GTPSA.mad_tpsa_poisbra!Method
mad_tpsa_poisbra!(a::Ptr{RTPSA}, b::Ptr{RTPSA}, c::Ptr{RTPSA}, nv::Cint)

Sets TPSA c to the poisson bracket of TPSAs a and b.

Input

  • a – Source TPSA a
  • b – Source TPSA b
  • nv – Number of variables in the TPSA

Output

  • c – Destination TPSA c
GTPSA.mad_tpsa_pow!Method
mad_tpsa_pow!(a::Ptr{RTPSA}, b::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets the destination TPSA c = a ^ b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a ^ b
GTPSA.mad_tpsa_powi!Method
mad_tpsa_powi!(a::Ptr{RTPSA}, n::Cint, c::Ptr{RTPSA})

Sets the destination TPSA c = a ^ n where n is an integer.

Input

  • a – Source TPSA a
  • n – Integer power

Output

  • c – Destination TPSA c = a ^ n
GTPSA.mad_tpsa_pown!Method
mad_tpsa_pown!(a::Ptr{RTPSA}, v::Cdouble, c::Ptr{RTPSA})

Sets the destination TPSA c = a ^ v where v is of double precision.

Input

  • a – Source TPSA a
  • v – "double" precision power

Output

  • c – Destination TPSA c = a ^ v
GTPSA.mad_tpsa_printMethod
mad_tpsa_print(t::Ptr{RTPSA}, name_::Cstring, eps_::Cdouble, nohdr_::Cint, stream_::Ptr{Cvoid})

Prints the TPSA coefficients with precision eps_. If nohdr_ is not zero, the header is not printed.

Input

  • t – TPSA to print
  • name_ – (Optional) Name of TPSA
  • eps_ – (Optional) Precision to output
  • nohdr_ – (Optional) If True, no header is printed
  • stream_ – (Optional) FILE pointer of output stream. Default is stdout
GTPSA.mad_tpsa_scanMethod
mad_tpsa_scan(stream_::Ptr{Cvoid})::Ptr{RTPSA}

Scans in a TPSA from the stream_.

Input

  • stream_ – (Optional) I/O stream from which to read the TPSA, default is stdin

Output

  • t – TPSA scanned from I/O stream_
GTPSA.mad_tpsa_scan_coef!Method
mad_tpsa_scan_coef!(t::Ptr{RTPSA}, stream_::Ptr{Cvoid})

Read TPSA coefficients into TPSA t. This should be used with mad_tpsa_scan_hdr for external languages using this library where the memory is managed NOT on the C side.

Input

  • stream_ – (Optional) I/O stream to read TPSA from, default is stdin

Output

  • t – TPSA with coefficients scanned from stream_
GTPSA.mad_tpsa_scan_hdrMethod
mad_tpsa_scan_hdr(kind_::Ref{Cint}, name_::Ptr{Cuchar}, stream_::Ptr{Cvoid})::Ptr{Desc}

Read TPSA header. Returns descriptor for TPSA given the header. This is useful for external languages using this library where the memory is managed NOT on the C side.

Input

  • kind_ – (Optional) Real or complex TPSA, or detect automatically if not provided.
  • name_ – (Optional) Name of TPSA
  • stream_ – (Optional) I/O stream to read TPSA from, default is stdin

Output

  • ret – Descriptor for the TPSA
GTPSA.mad_tpsa_scl!Method
mad_tpsa_scl!(a::Ptr{RTPSA}, v::Cdouble, c::Ptr{RTPSA})

Sets TPSA c to v*a.

Input

  • a – Source TPSA a
  • v – Scalar with double precision

Output

  • c – Destination TPSA c = v*a
GTPSA.mad_tpsa_sclord!Method
mad_tpsa_sclord!(t::Ptr{RTPSA}, r::Ptr{RTPSA}, inv::Cuchar, prm::Cuchar)

Scales all coefficients by order. If inv == 0, scales coefficients by order (derivation), else scales coefficients by 1/order (integration).

Input

  • t – Source TPSA
  • inv – Put order up, divide, scale by inv of value of order
  • prm – Parameters flag. If set to 0x0, the scaling excludes the order of the parameters in the monomials. Else, scaling is with total order of monomial

Output

  • r – Destination TPSA
GTPSA.mad_tpsa_set0!Method
mad_tpsa_set0!(t::Ptr{RTPSA}, a::Cdouble, b::Cdouble)

Sets the 0th order coefficient (scalar part of TPSA) according to coef[0] = a*coef[0] + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • a – Scaling of current 0th order value
  • b – Constant added to current 0th order value
GTPSA.mad_tpsa_seti!Method
mad_tpsa_seti!(t::Ptr{RTPSA}, i::Cint, a::Cdouble, b::Cdouble)

Sets the coefficient of monomial at index i to coef[i] = a*coef[i] + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • i – Index of monomial
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_tpsa_setm!Method
mad_tpsa_setm!(t::Ptr{RTPSA}, n::Cint, m::Vector{Cuchar}, a::Cdouble, b::Cdouble)

Sets the coefficient of monomial defined by byte array m to coef = a*coef + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as byte array
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_tpsa_setnam!Method
mad_tpsa_setnam!(t::Ptr{RTPSA}, nam::Cstring)

Sets the name of the RTPSA.

Input

  • t – TPSA
  • nam – Name to set for RTPSA
GTPSA.mad_tpsa_sets!Method
mad_tpsa_sets!(t::Ptr{RTPSA}, n::Cint, s::Cstring, a::Cdouble, b::Cdouble)

Sets the coefficient of monomial defined by string s to coef = a*coef + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • n – Length of monomial
  • s – Monomial as string
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_tpsa_setsm!Method
mad_tpsa_setsm!(t::Ptr{RTPSA}, n::Cint, m::Vector{Cint}, a::Cdouble, b::Cdouble)

Sets the coefficient of monomial defined by sparse monomial m to coef = a*coef + b. Does not modify other values in TPSA.

Input

  • t – TPSA
  • n – Length of monomial
  • m – Monomial as sparse monomial
  • a – Scaling of current coefficient
  • b – Constant added to current coefficient
GTPSA.mad_tpsa_setv!Method
mad_tpsa_setv!(t::Ptr{RTPSA}, i::Cint, n::Cint, v::Vector{Cdouble})

Vectorized setter of the coefficients for monomials with indices i..i+n. Useful for putting a matrix into a map.

Input

  • t – TPSA
  • i – Starting index of monomials to set coefficients
  • n – Number of monomials to set coefficients of starting at i
  • v – Array of coefficients for monomials i..i+n
GTPSA.mad_tpsa_setval!Method
mad_tpsa_setval!(t::Ptr{RTPSA}, v::Cdouble)

Sets the scalar part of the TPSA to v and all other values to 0 (sets the TPSA order to 0).

Input

  • t – TPSA to set to scalar
  • v – Scalar value to set TPSA
GTPSA.mad_tpsa_setvar!Method
mad_tpsa_setvar!(t::Ptr{RTPSA}, v::Cdouble, iv_::Cint, scl_::Cdouble)

Sets the 0th and 1st order values for the variables, and sets the rest of the variables to 0

Input

  • t – TPSA
  • v – 0th order value (coefficient)
  • iv_ – Variable index, optional if order of TPSA is 0 (behaves like mad_tpsa_setval then)
  • scl_ – 1st order variable value (typically will be 1)
GTPSA.mad_tpsa_sin!Method
mad_tpsa_sin!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the sin of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sin(a)
GTPSA.mad_tpsa_sinc!Method
mad_tpsa_sinc!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the sinc of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sinc(a)
GTPSA.mad_tpsa_sincos!Method
mad_tpsa_sincos!(a::Ptr{RTPSA}, s::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA s = sin(a) and TPSA c = cos(a)

Input

  • a – Source TPSA a

Output

  • s – Destination TPSA s = sin(a)
  • c – Destination TPSA c = cos(a)
GTPSA.mad_tpsa_sincosh!Method
mad_tpsa_sincosh!(a::Ptr{RTPSA}, s::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA s = sinh(a) and TPSA c = cosh(a)

Input

  • a – Source TPSA a

Output

  • s – Destination TPSA s = sinh(a)
  • c – Destination TPSA c = cosh(a)
GTPSA.mad_tpsa_sinh!Method
mad_tpsa_sinh!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the sinh of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sinh(a)
GTPSA.mad_tpsa_sinhc!Method
mad_tpsa_sinhc!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the sinhc of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sinhc(a)
GTPSA.mad_tpsa_sqrt!Method
mad_tpsa_sqrt!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the sqrt of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = sqrt(a)
GTPSA.mad_tpsa_sub!Method
mad_tpsa_sub!(a::Ptr{RTPSA}, b::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets the destination TPSA c = a - b

Input

  • a – Source TPSA a
  • b – Source TPSA b

Output

  • c – Destination TPSA c = a - b
GTPSA.mad_tpsa_tan!Method
mad_tpsa_tan!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the tan of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = tan(a)
GTPSA.mad_tpsa_tanh!Method
mad_tpsa_tanh!(a::Ptr{RTPSA}, c::Ptr{RTPSA})

Sets TPSA c to the tanh of TPSA a.

Input

  • a – Source TPSA a

Output

  • c – Destination TPSA c = tanh(a)
GTPSA.mad_tpsa_taylor!Method
mad_tpsa_taylor!(a::Ptr{RTPSA}, n::Cint, coef::Vector{Cdouble}, c::Ptr{RTPSA})

Computes the result of the Taylor series up to order n-1 with Taylor coefficients coef for the scalar value in a. That is, c = coef[0] + coef[1]*a_0 + coef[2]*a_0^2 + ... where a_0 is the scalar part of TPSA a.

Input

  • a – TPSA a
  • nOrder-1 of Taylor expansion, size of coef array
  • coef – Array of coefficients in Taylor s
  • c – Result
GTPSA.mad_tpsa_translate!Method
mad_tpsa_translate!(na::Cint, ma::Vector{Ptr{RTPSA}}, nb::Cint, tb::Vector{Cdouble}, mc::Vector{Ptr{RTPSA}})

Translates the expansion point of the map by the amount tb.

Input

  • na – Number of TPSAS in the map
  • ma – map ma
  • nb – Length of tb
  • tb – Vector of amount to translate for each variable

Output

  • mc – Map evaluated at the new point translated tb from the original evaluation point
GTPSA.mad_tpsa_uid!Method
mad_tpsa_uid!(t::Ptr{RTPSA}, uid_::Cint)::Cint

Sets the TPSA uid if uid_ != 0, and returns the current (previous if set) TPSA uid.

Input

  • t – TPSA
  • uid_uid to set in the TPSA if uid_ != 0

Output

  • ret – Current (previous if set) TPSA uid
GTPSA.mad_tpsa_unit!Method
mad_tpsa_unit!(x::Ptr{RTPSA}, r::Ptr{RTPSA})

Interpreting TPSA as a vector, gets the "unit vector", e.g. r = x/norm(x). May be useful for checking for convergence.

Input

  • x – Source TPSA x

Output

  • r – Destination TPSA r
GTPSA.mad_tpsa_vec2fld!Method
mad_tpsa_vec2fld!(na::Cint, a::Ptr{RTPSA}, mc::Vector{Ptr{RTPSA}})

Assuming the variables in the TPSA are canonically-conjugate, and ordered so that the canonically- conjugate variables are consecutive (q1, p1, q2, p2, ...), calculates the vector field (Hamilton's equations) from the passed Hamiltonian, defined as [da/dp1, -da/dq1, ...]

Input

  • na – Number of TPSA in mc consistent with number of variables in a
  • a – Hamiltonian as a TPSA

Output

  • mc – Vector field derived from a using Hamilton's equations
GTPSA.monoFunction
mono(v::Union{Integer, Vector{<:Union{<:Pair{<:Integer,<:Integer},<:Integer}}, Nothing}=nothing; param::Union{<:Integer,Nothing}=nothing, params::Union{Vector{<:Pair{<:Integer,<:Integer}}, Nothing}=nothing, use::Descriptor=GTPSA.desc_current)::TPS

Returns a TPS corresponding to a specific monomial, specified using the variable/parameter index, or monomial indexing-by-order OR monomial indexing-by-sparse monomial.

Input

  • v – An integer (for variable index), an array of orders for each variable (for indexing-by-order), or an array of pairs (sparse monomial)
  • param – (Keyword argument, optional) An integer for the parameter index
  • params – (Keyword argument, optional) An array of pairs for sparse-monomial indexing
  • use – (Keyword argument, optional) The descriptor to use to generate the monomial. Default is most recently-defined.

Examples: Variable/Parameter Index:

julia> d = Descriptor(3,10,2,10);

julia> mono(1)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1    0    0    0    0


julia> mono(2, use=d)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    1    0    0    0


julia> mono(param=2)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    0    0    0    1

Examples: Monomial Index-by-Order

julia> mono([1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1    0    0    0    0


julia> mono([0,1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    1    0    0    0


julia> mono([0,0,0,0,1], use=d)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    0    0    0    1


julia> mono([1,0,0,0,1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    2        1    0    0    0    1

Examples: Monomial Index-by-Sparse Monomial

julia> mono([1=>1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1    0    0    0    0


julia> mono([2=>1])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    1    0    0    0


julia> mono([1=>1], params=[2=>1], use=d)
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    2        1    0    0    0    1
GTPSA.parFunction
par(t::Union{TPS,ComplexTPS}, v::Union{Integer, Vector{<:Union{<:Pair{<:Integer,<:Integer},<:Integer, <:Any}}, Tuple{Vararg{Union{<:Integer,Pair{<:Integer,<:Integer},<:Colon}}}, Nothing}=nothing; param::Union{<:Integer,Nothing}=nothing, params::Union{Vector{<:Pair{<:Integer,<:Integer}}, Nothing}=nothing)

Extracts a polynomial from the TPS containing the specified monomial, and removes the monomial.

Input

  • v – An integer (for variable index), an array of orders for each variable (for indexing-by-order), or an array of pairs (sparse monomial)
  • param – (Keyword argument, optional) An integer for the parameter index
  • params – (Keyword argument, optional) An array of pairs for sparse-monomial indexing

Examples: Variable/Parameter Index:

julia> d = Descriptor(5, 10, 2, 10); x = vars(d); k = params(d);

julia> f = 2*x[1]^2*x[3] + 3*x[1]^2*x[2]*x[3]*x[4]^2*x[5]*k[1] + 6*x[3] + 5
TPS:
 Coefficient                Order   Exponent
  5.0000000000000000e+00      0      0   0   0   0   0   |   0   0
  6.0000000000000000e+00      1      0   0   1   0   0   |   0   0
  2.0000000000000000e+00      3      2   0   1   0   0   |   0   0
  3.0000000000000000e+00      8      2   1   1   2   1   |   1   0


julia> par(f, 3)
TPS:
 Coefficient                Order   Exponent
  6.0000000000000000e+00      0      0   0   0   0   0   |   0   0
  2.0000000000000000e+00      2      2   0   0   0   0   |   0   0
  3.0000000000000000e+00      7      2   1   0   2   1   |   1   0


julia> par(f, param=1)
TPS:
 Coefficient                Order   Exponent
  3.0000000000000000e+00      7      2   1   1   2   1   |   0   0

Examples: Monomial Index-by-Order

julia> par(f, [2,:,1])
TPS:
 Coefficient                Order   Exponent
  2.0000000000000000e+00      0      0   0   0   0   0   |   0   0
  3.0000000000000000e+00      5      0   1   0   2   1   |   1   0


julia> par(f, [2,0,1])
TPS:
 Coefficient                Order   Exponent
  2.0000000000000000e+00      0      0   0   0   0   0   |   0   0

Examples: Monomial Index-by-Sparse Monomial

julia> par(f, [1=>2, 3=>1])
TPS:
 Coefficient                Order   Exponent
  2.0000000000000000e+00      0      0   0   0   0   0   |   0   0
  3.0000000000000000e+00      5      0   1   0   2   1   |   1   0

  
julia> par(f, params=[1=>1])
TPS:
 Coefficient                Order   Exponent
  3.0000000000000000e+00      7      2   1   1   2   1   |   0   0
GTPSA.paramsFunction
params(use::Union{Descriptor,TPS,ComplexTPS}=GTPSA.desc_current)::Vector{TPS}

Returns TPSs corresponding to the parameters for the Descriptor specified by use. Default value is GTPSA.desc_current.

Input

  • use – (Optional) Specify which TPSA Descriptor to use, default is GTPSA.desc_current

Output

  • kVector containing unit TPSs corresponding to each parameter
GTPSA.pbMethod
pb(f::Union{TPS, ComplexTPS}, g::Union{TPS, ComplexTPS})

Assuming the variables in the TPSA are canonically-conjugate, and ordered so that the canonically- conjugate variables are consecutive (q₁, p₁, q₂, p₂, ...), computes the Poisson bracket of the scalar functions f and g. The Poisson bracket of two functions {f, g} is defined as Σᵢ (∂f/∂qᵢ)(∂g/∂pᵢ) - (∂g/∂qᵢ)(∂f/∂pᵢ).

Examples

julia> d = Descriptor(4,10);

julia> x = vars(d);

julia> f = (x[1]^2 + x[2]^2)/2 + (x[3]^2 + x[4]^2)/2;

julia> pb(f,x[1])
TPS:
  Coefficient              Order     Exponent
  -1.0000000000000000e+00    1        0    1    0    0


julia> pb(f,x[2])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        1    0    0    0


julia> pb(f,x[3])
TPS:
  Coefficient              Order     Exponent
  -1.0000000000000000e+00    1        0    0    0    1


julia> pb(f,x[4])
TPS:
  Coefficient              Order     Exponent
   1.0000000000000000e+00    1        0    0    1    0
GTPSA.ptinvMethod
ptinv(ma::Vector{<:Union{TPS,ComplexTPS}}, vars::Vector{<:Integer})

Partially-inverts the map ma, inverting only the variables specified by index in vars.

GTPSA.scalarMethod
scalar(t::Union{TPS,ComplexTPS})

Extracts the scalar part of the TPS. Equivalent to t[0] but this can be easily broadcasted.

GTPSA.translateMethod
translate(m::Vector{ComplexTPS}, x::Vector{<:Number})::Vector{ComplexTPS}

Translates the expansion point of the Vector of TPSs m by x.

GTPSA.translateMethod
translate(m::Vector{TPS}, x::Vector{<:Real})::Vector{TPS}

Translates the expansion point of the Vector of TPSs m by x.

GTPSA.varsFunction
vars(use::Union{Descriptor,TPS,ComplexTPS}=GTPSA.desc_current)::Vector{TPS}

Returns TPSs corresponding to the variables for the Descriptor specified by use. Default value is GTPSA.desc_current.

Input

  • use – (Optional) Specify which TPSA Descriptor to use, default is GTPSA.desc_current

Output

  • xVector containing unit TPSs corresponding to each variable
LinearAlgebra.normMethod
norm(t1::ComplexTPS)::Float64

Calculates the 1-norm of the ComplexTPS, which is the sum of the abs of all coefficients.

LinearAlgebra.normMethod
norm(t1::TPS)::Float64

Calculates the 1-norm of the TPS, which is the sum of the abs of all coefficients.

GTPSA.@FastGTPSAMacro
FastGTPSA(expr)

Macro speed up evaluation of mathematical expressions containing TPSs. @FastGTPSA is completely transparent to all other types, so it can be prepended to expressions while still maintaining type-generic code.

Example

julia> using GTPSA, BenchmarkTools

julia> d = Descriptor(3,5); x = vars(d);

julia> @btime $x[1]^3*sin($x[2])/log(2+$x[3])-exp($x[1]*$x[2])*im;
  2.114 μs (10 allocations: 160 bytes)

julia> @btime @FastGTPSA $x[1]^3*sin($x[2])/log(2+$x[3])-exp($x[1]*$x[2])*im;
  1.744 μs (1 allocation: 16 bytes)