Solver

The main symbolic solver for Symbolics.jl is symbolic_solve. Symbolic solving means that it only uses symbolic (algebraic) methods and outputs exact solutions.

Symbolics.symbolic_solveFunction
symbolic_solve(expr, x; dropmultiplicity=true, warns=true)

symbolic_solve is a function which attempts to solve input equations/expressions symbolically using various methods.

Arguments

  • expr: Could be a single univar expression in the form of a poly or multiple univar expressions or multiple multivar polys or a transcendental nonlinear function.

  • x: Could be a single variable or an array of variables which should be solved

  • dropmultiplicity (optional): Should the output be printed n times where n is the number of occurrence of the root? Say we have (x+1)^2, we then have 2 roots x = -1, by default the output is [-1], If dropmultiplicity is inputted as false, then the output is [-1, -1].

  • warns (optional): When invalid expressions or cases are inputted, should the solver warn you of such cases before returning nothing? if this is set to false, the solver returns nothing. By default, warns are set to true.

Supported input

The base solver (symbolic_solve) has multiple solvers which chooses from depending on the the type of input (multiple/uni var and multiple/single expression) only after ensuring that the input is valid.

The expressions inputted can contain parameters, which are assumed to be transcendental. A parameter "a" is transcendental if there exists no polynomial P with rational coefficients such that P(a) = 0. Check the examples section.

Currently, symbolic_solve supports

  • Linear and polynomial equations (with parameters)
  • Systems of linear and polynomials equations (without extra parameters, for now)
  • Equations with transcendental functions (with parameters)

Examples

solve_univar (uses factoring and analytic solutions up to degree 4)

Note

The package Nemo is needed in order to use solve_univar as well as solve_multipoly, so executing using Nemo as you will see in the following examples is necessary; otherwise, the function will throw an error.

julia> using Symbolics, Nemo;

julia> @variables x a b;

julia> expr = expand((x + b)*(x^2 + 2x + 1)*(x^2 - a))
-a*b - a*x - 2a*b*x - 2a*(x^2) + b*(x^2) + x^3 - a*b*(x^2) - a*(x^3) + 2b*(x^3) + 2(x^4) + b*(x^4) + x^5

julia> symbolic_solve(expr, x)
4-element Vector{Any}:
 -1
   -b
   (1//2)*√(4a)
   (-1//2)*√(4a)

julia> symbolic_solve(expr, x, dropmultiplicity=false)
5-element Vector{Any}:
 -1
 -1
   -b
   (1//2)*√(4a)
   (-1//2)*√(4a)
julia> symbolic_solve(x^2 + a*x + 6, x)
2-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 (1//2)*(-a + √(-24 + a^2))
 (1//2)*(-a - √(-24 + a^2))
julia> symbolic_solve(x^7 - 1, x)
2-element Vector{Any}:
  roots_of((1//1) + x + x^2 + x^3 + x^4 + x^5 + x^6, x)
 1

solve_multivar (uses Groebner basis and solve_univar to find roots)

Note

Similar to solve_univar, Groebner is needed for solve_multivar or to be fully functional.

julia> using Groebner

julia> @variables x y z
3-element Vector{Num}:
 x
 y
 z

julia> eqs = [x+y^2+z, z*x*y, z+3x+y]
3-element Vector{Num}:
 x + z + y^2
       x*y*z
  3x + y + z

julia> symbolic_solve(eqs, [x,y,z])
3-element Vector{Any}:
 Dict{Num, Any}(z => 0, y => 1//3, x => -1//9)
 Dict{Num, Any}(z => 0, y => 0, x => 0)
 Dict{Num, Any}(z => -1, y => 1, x => 0)
Note

If Nemo or Groebner are not imported when needed, the solver throws an error.

julia> using Symbolics

julia> @variables x y z;

julia> symbolic_solve(x+1, x)
ERROR: "Nemo is required. Execute `using Nemo` to enable this functionality."

julia> symbolic_solve([x+1, y], [x, y])
ERROR: "Groebner bases engine is required. Execute `using Groebner` to enable this functionality."

solve_multipoly (uses GCD between the input polys)

julia> symbolic_solve([x-1, x^3 - 1, x^2 - 1, (x-1)^20], x)
1-element Vector{BigInt}:
 1

ia_solve (solving by isolation and attraction)

julia> symbolic_solve(2^(x+1) + 5^(x+3), x)
1-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 (-slog(2) - log(complex(-1)) + 3slog(5)) / (slog(2) - slog(5))
julia> symbolic_solve(log(x+1)+log(x-1), x)
2-element Vector{SymbolicUtils.BasicSymbolic{BigFloat}}:
 (1//2)*√(8.0)
 (-1//2)*√(8.0)
julia> symbolic_solve(a*x^b + c, x)
((-c)^(1 / b)) / (a^(1 / b))

Evaluating output (converting to floats)

If you want to evaluate the exact expressions found by symbolic_solve, you can do the following:

julia> roots = symbolic_solve(2^(x+1) + 5^(x+3), x)
1-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 (-slog(2) - log(complex(-1)) + 3slog(5)) / (slog(2) - slog(5))

julia> Symbolics.symbolic_to_float.(roots)
1-element Vector{Complex{BigFloat}}:
 -4.512941594732059759689023145584186058252768936052415430071569066192919491762214 + 3.428598090438030380369414618548038962770087500755160535832807433942464545729382im
source

One other symbolic solver is symbolic_linear_solve which is limited compared to symbolic_solve as it only solves linear equations.

Symbolics.symbolic_linear_solveFunction
symbolic_linear_solve(eq, var; simplify, check) -> Any

Solve equation(s) eqs for a set of variables vars.

Assumes length(eqs) == length(vars)

Currently only works if all equations are linear. check if the expr is linear w.r.t vars.

Examples

julia> @variables x y
2-element Vector{Num}:
 x
 y

julia> Symbolics.symbolic_linear_solve(x + y ~ 0, x)
-y

julia> Symbolics.symbolic_linear_solve([x + y ~ 0, x - y ~ 2], [x, y])
2-element Vector{Float64}:
  1.0
 -1.0
source

symbolic_solve only supports symbolic, i.e. non-floating point computations, and thus prefers equations where the coefficients are integer, rational, or symbolic. Floating point coefficients are transformed into rational values and BigInt values are used internally with a potential performance loss, and thus it is recommended that this functionality is only used with floating point values if necessary. In contrast, symbolic_linear_solve directly handles floating point values using standard factorizations.

More technical details and examples

Technical details

The symbolic_solve function uses 4 hidden solvers in order to solve the user's input. Its base, solve_univar, uses analytic solutions up to polynomials of degree 4 and factoring as its method for solving univariate polynomials. The function's solve_multipoly uses GCD on the input polynomials then throws passes the result to solve_univar. The function's solve_multivar uses Groebner basis and a separating form in order to create linear equations in the input variables and a single high degree equation in the separating variable [1]. Each equation resulting from the basis is then passed to solve_univar. We can see that essentially, solve_univar is the building block of symbolic_solve. If the input is not a valid polynomial and can not be solved by the algorithm above, symbolic_solve passes it to ia_solve, which attempts solving by attraction and isolation [2]. This only works when the input is a single expression and the user wants the answer in terms of a single variable. Say log(x) - a == 0 gives us [e^a].

Symbolics.solve_univarFunction
solve_univar(expression, x; dropmultiplicity=true)

This solver uses analytic solutions up to degree 4 to solve univariate polynomials. It first handles the special case of the expression being of operation ^. E.g. math (x+2)^{20}. We solve this by removing the int 20, then solving the poly math x+2 on its own. If the parameter mult of the solver is set to true, we then repeat the found roots of math x+2 twenty times before returning the results to the user.

Step 2 is filtering the expression after handling this special case, and then factoring it using factor_use_nemo. We then solve all the factors outputted using the analytic methods implemented in the function get_roots and its children.

Arguments

  • expr: Single symbolics Num or SymbolicUtils.BasicSymbolic expression. This is equated to 0 and then solved. E.g. expr = x+2, we solve x+2 = 0

  • x: Single symbolics variable

  • dropmultiplicity (optional): Print repeated roots or not?

  • strict (optional): Bool that enables/disables strict assert if input expression is a univariate polynomial or not. If strict=true and expression is not a polynomial, solve_univar throws an assertion error.

Examples

source
Missing docstring.

Missing docstring for Symbolics.solve_multivar. Check Documenter's build log for details.

Symbolics.ia_solveFunction
ia_solve(lhs, var; kwargs...)

This function attempts to solve transcendental functions by first checking the "smart" number of occurrences in the input LHS. By smart here we mean that polynomials are counted as 1 occurrence. for example x^2 + 2x is 1 occurrence of x. So we abstract all occurrences of x's as polynomials. Say: log(x+1) + x^2 is seen as log(f(x)) + g(x) so there are 2 occurrences of x. If there is only 1 occurrence of x in an input expression, isolate is called.

Isolate reverses all operations applied on the occurrence of x until we have f(x) = some constant then we can solve this using our polynomial solvers.

If more than 1 occurrence of x is found, ia_solve attempts to attract the occurrences of x in order to reduce these occurrences to 1. For example, log(x+1) + log(x-1) can be converted to log(x^2 - 1) which now could be isolated using Isolate.

attract(lhs, var) currently uses 4 techniques for attraction.

  • Log addition: log(f(x)) + log(g(x)) => log(h(x))
  • Exponential simplification: a*b^(f(x)) + c*d^(g(x)) => f(x) * log(b) - g(x) * log(d) + log(-a/c). And now this is actually 1 occurrence of x since f(x) and g(x) are just multiplied by constants not wrapped in some operation.
  • Trig simplification: this bruteforces multiple trig identities and doesn't detect them before hand.
  • Polynomialization: as a last resort, attract attempts to polynomialize the expression. Say sin(x+2)^2 + sin(x+2) + 10 is converted to X^2 + X + 10, we then solve this using our polynomial solver, and afterwards, isolate sin(x+2) = the roots found by solve for X^2 + X + 10

After attraction, we check the number of occurrences again, and if its 1, we isolate, if not, we throw an error to tell the user that this is currently unsolvable by our covered techniques.

Arguments

  • lhs: a Num/SymbolicUtils.BasicSymbolic
  • var: variable to solve for.

Keyword arguments

  • warns = true: Whether to emit warnings for unsolvable expressions.
  • complex_roots = true: Whether to consider complex roots of x ^ n ~ y, where n is an integer.
  • periodic_roots = true: If true, isolate f(x) ~ y as x ~ finv(y) + n * period where is_periodic(f) == true, finv = left_inverse(f) and period = fundamental_period(f). n is a new anonymous symbolic variable.

Examples

julia> solve(a*x^b + c, x)
((-c)^(1 / b)) / (a^(1 / b))
julia> solve(2^(x+1) + 5^(x+3), x)
1-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 (-log(2) + 3log(5) - log(complex(-1))) / (log(2) - log(5))
julia> solve(log(x+1)+log(x-1), x)
2-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 (1//2)*RootFinding.ssqrt(8.0)
 (-1//2)*RootFinding.ssqrt(8.0)
julia> expr = sin(x+2)^2 + sin(x+2) + 10
10 + sin(2 + x) + sin(2 + x)^2

julia> RootFinding.ia_solve(expr, x)
[ Info: var"##230" ϵ Ζ: e.g. 0, 1, 2...
[ Info: var"##234" ϵ Ζ: e.g. 0, 1, 2...
2-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 -2 + π*2var"##230" + asin((1//2)*(-1 + RootFinding.ssqrt(-39)))
 -2 + π*2var"##234" + asin((1//2)*(-1 - RootFinding.ssqrt(-39)))

All transcendental functions for which left_inverse is defined are supported. To enable ia_solve to handle custom transcendental functions, define an inverse or left inverse. If the function is periodic, is_periodic and fundamental_period must be defined. If the function imposes certain conditions on its input or output (for example, log requires that its input be positive) define ia_conditions!.

See also: left_inverse, inverse, is_periodic, fundamental_period, ia_conditions!.

References

source
Symbolics.ia_conditions!Function
ia_conditions!(f, lhs, rhs::Vector{Any}, conditions::Vector{Tuple})

If f is a left-invertible function, lhs and rhs[i] are univariate functions and f(lhs) ~ rhs[i] for all i in eachindex(rhss), push to conditions all the relevant conditions on lhs or rhs[i]. Each condition is of the form (sym, op) where sym is an expression involving lhs and/or rhs[i] and op is a binary relational operator. The condition op(sym, 0) is then required to be true for the equation f(lhs) ~ rhs[i] to be valid.

For example, if f = log, lhs = x and rhss = [y, z] then the condition x > 0 must be true. Thus, (lhs, >) is pushed to conditions. Similarly, if f = sqrt, rhs[i] >= 0 must be true for all i, and so (y, >=) and (z, >=) will be appended to conditions.

source
Symbolics.is_periodicFunction
is_periodic(f)

Return true if f is a single-input single-output periodic function. Return false by default. If is_periodic(f) == true, then fundamental_period(f) must also be defined.

See also: fundamental_period

source

Nice examples

using Symbolics, Nemo;
@variables x;
Symbolics.symbolic_solve(9^x + 3^x ~ 8, x)
2-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 slog(-(1//2) + (1//2)*√(33)) / slog(3)
 slog(-(1//2) - (1//2)*√(33)) / slog(3)
@variables x y z;
Symbolics.symbolic_linear_solve(2//1*x + y - 2//1*z ~ 9//1*x, 1//1*x)

\[ \begin{equation} \frac{1}{7} \left( y - 2 z \right) \end{equation} \]

using Groebner;
@variables x y z;

eqs = [x^2 + y + z - 1, x + y^2 + z - 1, x + y + z^2 - 1]
Symbolics.symbolic_solve(eqs, [x,y,z])
5-element Vector{Any}:
 Dict{Num, Any}(z => 1, y => 0, x => 0)
 Dict{Num, Any}(z => -1 + √(2), y => -1 + √(2), x => -1 + √(2))
 Dict{Num, Any}(z => -1 - √(2), y => -1 - √(2), x => -1 - √(2))
 Dict{Num, Any}(z => 0, y => 0, x => 1)
 Dict{Num, Any}(z => 0, y => 1, x => 0)

Feature completeness

  • [x] Linear and polynomial equations
  • [x] Systems of linear and polynomial equations
  • [x] Some transcendental functions
  • [x] Systems of linear equations with parameters (via symbolic_linear_solve)
  • [ ] Equations with radicals
  • [x] Systems of polynomial equations with parameters and positive dimensional systems
  • [ ] Inequalities

Expressions we can not solve (but aim to)

# Mathematica

In[1]:= Reduce[x^2 - x - 6 > 0, x]
Out[1]= x < -2 || x > 3

In[2]:= Reduce[x+a > 0, x]
Out[2]= a \[Element] Reals && x > -a

In[3]:= Solve[x^(x)  + 3 == 0, x]
Out[3]= {{x -> (I \[Pi] + Log[3])/ProductLog[I \[Pi] + Log[3]]}}

References