ModelKit
ModelKit is the symbolic input and modeling language of HomotopyContinuation.jl. It is designed such that you can easily create an efficient formulation of your problem.
Expressions and Variables
HomotopyContinuation.ModelKit.Expression — TypeExpression <: NumberA symbolic expression.
julia> expr = (Variable(:x) + 1)^2
(1 + x)^2
julia> Expression(2)
2
julia> Expression(Variable(:x))
xHomotopyContinuation.ModelKit.Variable — TypeVariable(name::Union{String,Symbol}, indices...) <: NumberA data structure representing a variable.
julia> Variable(:a)
a
julia> Variable(:x, 1)
x₁
julia> Variable(:x, 10, 5)
x₁₀₋₅Equality and ordering
Variables are identified by their name and indices. That is, two variables are equal if and only if they have the same name and indices.
julia> Variable(:a) == Variable(:a)
true
julia> Variable(:a, 1) == Variable(:a, 2)
falseSimilarly, variables are first ordered lexicographically by their name and then by their indices.
julia> Variable(:a, 1) < Variable(:a, 2)
true
julia> Variable(:a, 1) < Variable(:b, 1)
true
julia> a = [Variable(:a, i, j) for i in 1:2, j in 1:2]
2×2 Array{Variable,2}:
a₁₋₁ a₁₋₂
a₂₋₁ a₂₋₂
julia> sort(vec(a))
4-element Array{Variable,1}:
a₁₋₁
a₂₋₁
a₁₋₂
a₂₋₂HomotopyContinuation.ModelKit.@var — Macro@var variable1 variable2 ...Declare variables with the given names and automatically create the variable bindings. The macro supports indexing notation to create Arrays of variables.
Examples
julia> @var a b x[1:2] y[1:2,1:3]
(a, b, Variable[x₁, x₂], Variable[y₁₋₁ y₁₋₂ y₁₋₃; y₂₋₁ y₂₋₂ y₂₋₃])
julia> a
a
julia> b
b
julia> x
2-element Array{Variable,1}:
x₁
x₂
julia> y
2×3 Array{Variable,2}:
y₁₋₁ y₁₋₂ y₁₋₃
y₂₋₁ y₂₋₂ y₂₋₃HomotopyContinuation.ModelKit.@unique_var — Macro@unique_var variable1 variable2This is similar to @var with the only difference that the macro automatically changes the names of the variables to ensure uniqueness. However, the binding is still to the declared name. This is useful to ensure that there are no name collisions.
Examples
julia> @unique_var a b
(a#591, b#592)
julia> a
a#591
julia> b
b#592MultivariatePolynomials.variables — Methodvariables(prefix::Union{Symbol,String}, indices...)Create an Array of variables with the given prefix and indices. The expression @var x[1:3, 1:2] is equivalent to x = variables(:x, 1:3, 1:2).
Methods
MultivariatePolynomials.coefficients — Methodcoefficients(f::Expression, vars::AbstractVector{Variable}; expanded = false)Return all coefficients of the given polynomial f for the given variables vars. If expanded = true then this assumes that the expression f is already expanded, e.g., with expand.
HomotopyContinuation.ModelKit.coeffs_as_dense_poly — Functioncoeffs_as_dense_poly(f, vars, d; homogeneous = false)Given a polynomial f this returns a vector c of coefficients such that subs(dense_poly(vars, d; homogeneous = homogeneous), c) == f.
Example
@var x[1:3]
f, c = dense_poly(x, 3, coeff_name = :c)
g = x[1]^3+x[2]^3+x[3]^3-1
gc = coeffs_as_dense_poly(g, x, 3)
subs(f, c => gc) == gtrueMultivariatePolynomials.degree — Methoddegree(f::Expression, vars = variables(f); expanded = false)Compute the degree of the expression f in vars. Unless expanded is true the expression is first expanded.
HomotopyContinuation.ModelKit.degrees — Methoddegrees(f::AbstractVector{Expression}, vars = variables(f); expanded = false)Compute the degrees of the expressions f in vars. Unless expanded is true the expressions are first expanded.
MultivariatePolynomials.differentiate — Methoddifferentiate(expr::Expression, var::Variable, k = 1)
differentiate(expr::AbstractVector{Expression}, var::Variable, k = 1)Compute the k-th derivative of expr with respect to the given variable var.
differentiate(expr::Expression, vars::AbstractVector{Variable})Compute the partial derivatives of expr with respect to the given variable variables vars. Retuns a Vector containing the partial derivatives.
differentiate(exprs::AbstractVector{Expression}, vars::AbstractVector{Variable})Compute the partial derivatives of exprs with respect to the given variable variables vars. Returns a Matrix where the each row contains the partial derivatives for a given expression.
Base.conj — Functionconj(A::AbstractArray)Return an array containing the complex conjugate of each entry in array A.
Equivalent to conj.(A), except that when eltype(A) <: Real A is returned without copying, and that when A has zero dimensions, a 0-dimensional array is returned (rather than a scalar).
Examples
julia> conj([1, 2im, 3 + 4im])
3-element Vector{Complex{Int64}}:
1 + 0im
0 - 2im
3 - 4im
julia> conj(fill(2 - im))
0-dimensional Array{Complex{Int64}, 0}:
2 + 1imconj(z)Compute the complex conjugate of a complex number z.
Examples
julia> conj(1 + 3im)
1 - 3imconj(x::AbstractPolynomialLike)Return the complex conjugate of x by applying conjugation to all coefficients and variables.
conj(x::AbstractVector{<:AbstractMonomial})Return the complex conjugate of x by applying conjugation to monomials.
conj(::AbstractPolynomial)Returns the complex conjugate of the polynomial
conj(P::IntegerPartition) returns the Ferrer's conjugate of P. Also available as P'.
Base.conj(expr::Expression)This command returns back expr. This is intended, because complex conjugation is not algebraic.
Example
julia> @var x y ;
julia> f = x + im*y
julia> conj(f)
x + im*yPolynomials with conjugated coefficients.
If f is a polynomial, the polynomial obtained from f by conjugating its coefficients can be obtained as in the following example.
julia> @var x y ;
julia> f = x + im*y
julia> vars = variables(f)
julia> e, c = exponents_coefficients(f, vars)
julia> poly_from_exponents_coefficients(e, conj.(c), vars)
x - im*yHomotopyContinuation.ModelKit.dense_poly — Functiondense_poly(vars::AbstractVector{Variable}, d::Integer;
homogeneous = false,
coeff_name::Symbol = gensym(:c))Create a dense polynomial of degree d in the given variables variables where each coefficient is a parameter. Returns a tuple with the first argument being the polynomial and the second the parameters.
julia> @var x y;
julia> f, c = dense_poly([x, y], 2, coeff_name = :q);
julia> f
q₆ + x*q₄ + x^2*q₁ + y*q₅ + y^2*q₃ + x*y*q₂
julia> c
6-element Array{Variable,1}:
q₁
q₂
q₃
q₄
q₅
q₆HomotopyContinuation.ModelKit.evaluate — Methodevaluate(expr::Expression, subs...)
evaluate(expr::AbstractArray{Expression}, subs...)Evaluate the given expression.
Example
julia> @var x y;
julia> evaluate(x^2, x => 2)
4
julia> evaluate(x * y, [x,y] => [2, 3])
6
julia> evaluate([x^2, x * y], [x,y] => [2, 3])
2-element Array{Int64,1}:
4
6
# You can also use the callable syntax
julia> [x^2, x * y]([x,y] => [2, 3])
2-element Array{Int64,1}:
4
6HomotopyContinuation.ModelKit.expand — Functionexpand(e::Expression)Expand a given expression.
julia> @var x y
(x, y)
julia> expand((x + y) ^ 2)
2*x*y + x^2 + y^2HomotopyContinuation.ModelKit.exponents_coefficients — Functionexponents_coefficients(f::Expression, vars::AbstractVector{Variable}; expanded = false)Return a matrix M containing the exponents for all occuring terms (one term per column) and a vector c containing the corresponding coefficients. Expands the given expression f unless expanded = true. Throws a PolynomialError if a rational expression is encountered.
Missing docstring for poly_from_exponents_coefficients. Check Documenter's build log for details.
HomotopyContinuation.ModelKit.horner — Functionhorner(f::Expression, vars = variables(f))Rewrite f using a multi-variate horner schema.
Example
julia> @var u v c[1:3]
(u, v, Variable[c₁, c₂, c₃])
julia> f = c[1] + c[2] * v + c[3] * u^2 * v^2 + c[3]u^3 * v
c₁ + v*c₂ + u^2*v^2*c₃ + u^3*v*c₃
julia> horner(f)
c₁ + v*(c₂ + u^3*c₃ + u^2*v*c₃)MultivariatePolynomials.nvariables — Methodnvariables(expr::Expression; parameters = Variable[])
nvariables(exprs::AbstractVector{Expression}; parameters = Variable[])Obtain the number of variables used in the given expression not counting the the ones declared in parameters.
MultivariatePolynomials.monomials — Methodmonomials(variables::AbstractVector, d::Integer; affine = true)Create all monomials of a given degree in the given variables.
julia> @var x y
(x, y)
julia> monomials([x,y], 2; affine = false)
3-element Array{Operation,1}:
x ^ 2
x * y
y ^ 2MultivariatePolynomials.subs — Methodsubs(expr::Expression, subsitutions::Pair...)
subs(exprs::AbstractVector{<:Expression}, subsitutions::Pair...)Apply the given substitutions to the given expressions.
Examples
@var x y
julia> subs(x^2, x => y)
y ^ 2
julia> subs(x * y, [x,y] => [x+2,y+2])
(x + 2) * (y + 2)
julia> subs([x + y, x^2], x => y + 2, y => x + 2)
2-element Array{Expression,1}:
4 + x + y
(2 + y)^2
# You can also use the callable syntax
julia> (x * y)([x,y] => [x+2,y+2])
(x + 2) * (y + 2)HomotopyContinuation.ModelKit.rand_poly — Functionrand_poly(T = ComplexF64, vars::AbstractVector{Variable}, d::Integer; homogeneous::Bool = false)Create a random dense polynomial of degree d in the given variables variables. Each coefficient is sampled independently via randn(T).
julia> @var x y;
julia> rand_poly(Float64, [x, y], 2)
0.788764085756728 - 0.534507647623108*x - 0.778441366874946*y -
0.128891763280247*x*y + 0.878962738754971*x^2 + 0.550480741774464*y^2HomotopyContinuation.ModelKit.to_dict — Functionto_dict(expr::Expression, vars::AbstractVector{Variable})Return the coefficients of expr w.r.t. the given variables vars. Assumes that expr is expanded and representing a polynomial. Throws a PolynomialError if a rational expression is encountered.
HomotopyContinuation.ModelKit.to_number — Functionto_number(x::Expression)Tries to unpack the Expression x to a native number type.
```julia-repl julia> x = to_number(Expression(2)) 2
julia> typeof(x) Int64
MultivariatePolynomials.variables — Methodvariables(expr::Expression; parameters = Variable[])
variables(exprs::AbstractVector{Expression}; parameters = Variable[])Obtain all variables used in the given expression up to the ones declared in parameters.
Example
julia> @var x y a;
julia> variables(x^2 + y)
2-element Array{Variable,1}:
x
y
julia> variables([x^2 + a, y]; parameters = [a])
2-element Array{Variable,1}:
x
ySystem
HomotopyContinuation.ModelKit.System — TypeSystem(exprs::AbstractVector{Expression};
variables = variables(exprssion),
parameters = Variable[])Create a system from the given Expressions exprs. The variables determine also the variable ordering. The parameters argument allows to declare certain Variables as parameters.
System(support::AbstractVector{<:AbstractMatrix{<:Integer}},
coefficients::AbstractVector{<:AbstractVector};
variables,
parameters = Variable[])Create a system from the given support and coefficients.
Examples
julia> @var x y;
julia> F = System([x^2, y^2]; variables = [y, x])
System of length 2
2 variables: y, x
x^2
y^2
# Systems are callable.
# This evaluates F at y=2 and x=3
julia> F([2, 3])
2-element Array{Int64,1}:
9
4It is also possible to declare parameters.
julia> @var x y a b;
julia> F = System([x^2 + a, y^2 + b]; variables = [y, x], parameters = [a, b])
System of length 2
2 variables: y, x
2 parameters: a, b
a + x^2
b + y^2
julia> F([2, 3], [5, -2])
2-element Array{Int64,1}:
14
2Missing docstring for evaluate(F::System, x, p = nothing). Check Documenter's build log for details.
HomotopyContinuation.ModelKit.jacobian — Methodjacobian(F::System)Computes the symbolic Jacobian of the system F.
Example
@var x y
F = System([x^2 + 3y, (y*x+1)^3])
jacobian(F)2×2 Array{Expression,2}:
2*x 3
3*y*(1 + x*y)^2 3*x*(1 + x*y)^2HomotopyContinuation.ModelKit.jacobian — Functionjacobian(F::System, x, p = nothing)Evaluates the Jacobian of the system F at (x, p).
Example
@var x y
F = System([x^2 + 3y, (y*x+1)^3])
jacobian(F, [2, 3])2×2 Array{Int32,2}:
4 3
441 294HomotopyContinuation.ModelKit.degrees — Methoddegrees(F::System)Return the degrees of the given system.
HomotopyContinuation.ModelKit.expressions — Methodexpressions(F::System)Returns the expressions of the given system F.
HomotopyContinuation.ModelKit.optimize — Methodoptimize(F::System)Optimize the evaluation cost of the given system F. This applies a multivariate horner schema to the given expressions. See also horner.
Example
julia> f
System of length 4
4 variables: z₁, z₂, z₃, z₄
z₁ + z₂ + z₃ + z₄
z₂*z₁ + z₂*z₃ + z₃*z₄ + z₄*z₁
z₂*z₃*z₁ + z₂*z₃*z₄ + z₂*z₄*z₁ + z₃*z₄*z₁
-1 + z₂*z₃*z₄*z₁
julia> optimize(f)
System of length 4
4 variables: z₁, z₂, z₃, z₄
z₁ + z₂ + z₃ + z₄
(z₂ + z₄)*z₁ + (z₂ + z₄)*z₃
z₁*(z₃*z₄ + (z₃ + z₄)*z₂) + z₂*z₃*z₄
-1 + z₂*z₃*z₄*z₁HomotopyContinuation.ModelKit.multi_degrees — Methodmulti_degrees(F::System)Return the degrees with respect to the given variable groups.
HomotopyContinuation.ModelKit.nparameters — Methodnparameters(F::System)Returns the number of parameters of the given system F.
MultivariatePolynomials.nvariables — Methodnvariables(F::System)Returns the number of variables of the given system F.
HomotopyContinuation.ModelKit.parameters — Methodparameters(F::System)Returns the parameters of the given system F.
HomotopyContinuation.ModelKit.support_coefficients — Methodsupport_coefficients(F::System)Return the support of the system and the corresponding coefficients.
MultivariatePolynomials.variables — Methodvariables(F::System)Returns the variables of the given system F.
HomotopyContinuation.ModelKit.variable_groups — Methodvariable_groups(F::System)Returns the variable groups of the given system F.
Homotopy
HomotopyContinuation.ModelKit.Homotopy — TypeHomotopy(exprs, vars, t, parameters = Variable[])Create a homotopy H(vars,t) from the given Expressions exprs where vars are the given variables and t is the dedicated variable parameterizing the family of systems. The parameters argument allows to declare certain Variables as parameters.
Example
julia> @var x y t;
julia> H = Homotopy([x + t, y + 2t], [y, x], t)
Homotopy in t of length 2
2 variables: y, x
t + x
2*t + y
julia> H([2, 3], 0)
2-element Array{Int64,1}:
3
2
julia> H([2, 3], 1)
2-element Array{Int64,1}:
4
4It is also possible to declare additional variables.
julia> @var x y t a b;
julia> H = Homotopy([x^2 + t*a, y^2 + t*b], [x, y], t, [a, b])
Homotopy in t of length 2
2 variables: x, y
2 parameters: a, b
a*t + x^2
b*t + y^2
julia> H([2, 3], 1, [5, 2])
2-element Array{Int64,1}:
9
11HomotopyContinuation.ModelKit.expressions — Methodexpressions(H::Homotopy)Returns the expressions of the given homotopy H.
HomotopyContinuation.ModelKit.nparameters — Methodnparameters(H::Homotopy)Returns the number of parameters of the given homotopy H.
MultivariatePolynomials.nvariables — Methodnvariables(H::Homotopy)Returns the number of variables of the given homotopy H.
HomotopyContinuation.ModelKit.parameters — Methodparameters(H::Homotopy)Returns the parameters of the given homotopy H.
MultivariatePolynomials.variables — Methodvariables(H::Homotopy)Returns the variables of the given homotopy H.