LACE.jl functions
The following functions are wrapped from the LACE.jl
package, which interfaces with the C-library LACE
. These functions transform the input polynomials and points into C-compatible objects, perform computations using interval arithmetic provided by MPFI
, and return results wrapped in Julia objects using MPFI.jl
as BigInterval
.
LACE.interval_classic_newton
— Functioninterval_classic_newton(sys::Vector{<:Union{AbstractAlgebra.Generic.MPoly, QQMPolyRingElem, ZZMPolyRingElem}},
initial_point::Vector{U},
input_precision::Int64=63,
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
verbose::Cint=Cint(1)) where {U<:Union{BigFloat,BigInterval}, T}
Performs the classical interval Newton method for solving a system of equations using AbstractAlgebra or Nemo polynomials.
Description
This function implements the classical interval Newton method through the following steps:
Input Processing:
- Converts the polynomial system to a C-object with interval coefficients at specified input precision
- Copies the initial point to C while preserving its original precision (BigFloat or BigInterval)
Solution Process:
- Performs the Kantorovich test to verify conditions for solution existence and uniqueness
- Executes the interval Newton method with the specified working precision
- Iteratively refines the solution until target precision is reached or max iterations exceeded
Output:
- Returns a tuple containing:
- Status code indicating success or failure
- Solution intervals rounded to output precision
- Kantorovich test constants
- Returns a tuple containing:
Arguments
sys::Vector
: A vector of polynomials from AbstractAlgebra or Nemo representing the system of equationsinitial_point::Vector{U}
: Starting point for the Newton method, where U is either BigFloat or BigIntervalinput_precision::Int64
: Precision (in bits) for input polynomial coefficients. Default is 63work_precision::Int64
: Precision for intermediate calculations. Default is 63target_precision::Int64
: Desired solution precision (interval length of 2^(-target_precision)). Default is 2output_precision::Int64
: Precision for rounding final result. Default is 2max_nb_loops::Int64
: Maximum number of iterations. Default is 10verbose::Cint
: Verbosity level for logging. Default is 1
Returns
Tuple (status, solution, constants)
where:
status::Int
: Result code (0: singular Jacobian, 1: solution found, 2: max precision reached, 3: max iterations reached)solution::Vector{BigInterval}
: Computed solution intervalsconstants::Tuple
: Kantorovich test constants
interval_classic_newton(f::Union{AbstractAlgebra.Generic.MPoly, QQMPolyRingElem, ZZMPolyRingElem},
initial_point::U,
input_precision::Int64=63,
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
verbose::Cint=Cint(1)) where {U<:Union{BigFloat,BigInterval},T}
Single equation version of interval Newton method for AbstractAlgebra or Nemo polynomials.
Arguments
f
: A single polynomial from AbstractAlgebra or Nemoinitial_point::U
: Starting point (BigFloat or BigInterval)- Other arguments as above
Returns
Same as above, but for a single equation
interval_classic_newton(sys::Vector{DynamicPolynomials.Polynomial{T, U, BigInterval}},
initial_point::Vector{V},
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
verbose::Cint=Cint(1)) where {V<:Union{BigFloat,BigInterval}, T, U}
Performs the classical interval Newton method for a system of DynamicPolynomials.
Description
This function implements the same classical interval Newton method as above, with one key difference:
- The system is already in interval form (DynamicPolynomials with BigInterval coefficients)
- No input precision parameter is needed as the system is copied with the maximum precision of its coefficients
Arguments
sys::Vector
: A vector of DynamicPolynomials with BigInterval coefficientsinitial_point::Vector{V}
: Starting point (BigFloat or BigInterval)work_precision::Int64
: Precision for intermediate calculations. Default is 63target_precision::Int64
: Desired solution precision. Default is 2output_precision::Int64
: Precision for rounding final result. Default is 2max_nb_loops::Int64
: Maximum number of iterations. Default is 10verbose::Cint
: Verbosity level. Default is 1
Returns
Same as above
interval_classic_newton(f::DynamicPolynomials.Polynomial{T, U, BigInterval},
initial_point::V,
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
verbose::Cint=Cint(1)) where {V<:Union{BigFloat,BigInterval}, T, U}
Single equation version of interval Newton method for DynamicPolynomials.
Arguments
f
: A single DynamicPolynomial with BigInterval coefficientsinitial_point::V
: Starting point (BigFloat or BigInterval)- Other arguments as above
Returns
Same as above, but for a single equation
interval_classic_newton(sys_ptr::Ptr{Cvoid},
initial_point::Vector{T},
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
verbose::Cint=Cint(1)) where {T<:Union{BigFloat,BigInterval}}
Low-level version of the classical interval Newton method that works directly with C pointers.
Arguments
sys_ptr::Ptr{Cvoid}
: A C pointer to the polynomial system (typically obtained fromsys_get_C_ptr
)initial_point::Vector{T}
: Starting point (BigFloat or BigInterval)work_precision::Int64
: Precision for intermediate calculations. Default is 63target_precision::Int64
: Desired solution precision. Default is 2output_precision::Int64
: Precision for rounding final result. Default is 2max_nb_loops::Int64
: Maximum number of iterations. Default is 10verbose::Cint
: Verbosity level. Default is 1
Returns
Tuple (status, solution, constants)
where:
status::Int
: Result code (0: singular Jacobian, 1: solution found, 2: max precision reached, 3: max iterations reached)solution::Vector{BigInterval}
: Computed solution intervalsconstants::Tuple
: Kantorovich test constants
Notes
- This is an internal function that should typically be called through the higher-level interfaces
- The system pointer must be properly initialized using
sys_get_C_ptr
Examples
# Using AbstractAlgebra
using AbstractAlgebra
R, (x,y) = polynomial_ring(QQ, ["x", "y"])
f = x^2 + y^2 - 1
g = x - y
sys = [f, g]
initial = [BigFloat(0.5), BigFloat(0.5)]
status, sol, consts = LACE.interval_classic_newton(sys, initial)
# Using DynamicPolynomials
using DynamicPolynomials, MPFI
@polyvar x y
f = BigInterval(1)*x^2 + BigInterval(1)*y^2 - BigInterval(1.0, 1.1)
g = x - y
sys = [f, g]
initial = [BigInterval(0.5), BigInterval(0.5)]
status, sol, consts = LACE.interval_classic_newton(sys, initial)
LACE.interval_newton
— Functioninterval_newton(sys::Vector{<:Union{AbstractAlgebra.Generic.MPoly, QQMPolyRingElem, ZZMPolyRingElem}},
initial_point::Vector{W},
input_precision::Int64=63,
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
refine::Int32=Int32(1),
verbose::Cint=Cint(1)) where {W<:Union{BigFloat,BigInterval}, T}
Performs the interval Newton method for solving a system of equations using AbstractAlgebra or Nemo polynomials. The refinement strategy (controlled by refine
) affects how the method handles successive iterations:
- With
refine=1
, the method always takes intersections, which is appropriate when the initial interval is known to contain the solution - With
refine=0
, the method uses a heuristic approach that only performs intersection when the new solution interval contains the center of the previous interval.
Description
This function implements the standard interval Newton method through the following steps:
Input Processing:
- Converts the polynomial system to a C-object with interval coefficients at specified input precision
- Copies the initial point to C while preserving its original precision (BigFloat or BigInterval)
Solution Process:
- Performs the Kantorovich test to verify conditions for solution existence and uniqueness
- Executes the interval Newton method with intersection operations between iterations
- Iteratively refines the solution until target precision is reached or max iterations exceeded
Output:
- Returns a tuple containing:
- Status code indicating success or failure
- Solution intervals rounded to output precision
- Kantorovich test constants
- Returns a tuple containing:
Arguments
sys::Vector
: A vector of polynomials from AbstractAlgebra or Nemo representing the system of equationsinitial_point::Vector{W}
: Starting point for the Newton method, where W is either BigFloat or BigIntervalinput_precision::Int64
: Precision (in bits) for input polynomial coefficients. Default is 63work_precision::Int64
: Precision for intermediate calculations. Default is 63target_precision::Int64
: Desired solution precision (interval length of 2^(-target_precision)). Default is 2output_precision::Int64
: Precision for rounding final result. Default is 2max_nb_loops::Int64
: Maximum number of iterations. Default is 10refine::Int32
: Flag to control refinement behavior. Default is 1verbose::Cint
: Verbosity level for logging. Default is 1
Returns
Tuple (status, solution, constants)
where:
status::Int
: Result code (0: singular Jacobian, 1: solution found, 2: max precision reached, 3: max iterations reached)solution::Vector{BigInterval}
: Computed solution intervalsconstants::Tuple
: Kantorovich test constants
interval_newton(f::Union{AbstractAlgebra.Generic.MPoly, QQMPolyRingElem, ZZMPolyRingElem},
initial_point::W,
input_precision::Int64=63,
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
refine::Int32=Int32(1),
verbose::Cint=Cint(1)) where {W<:Union{BigFloat,BigInterval},T}
Single equation version of interval Newton method for AbstractAlgebra or Nemo polynomials.
Arguments
f
: A single polynomial from AbstractAlgebra or Nemoinitial_point::W
: Starting point (BigFloat or BigInterval)- Other arguments as above
Returns
Same as above, but for a single equation
interval_newton(sys::Vector{DynamicPolynomials.Polynomial{T, U, BigInterval}},
initial_point::Vector{V},
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
refine::Int32=Int32(1),
verbose::Cint=Cint(1)) where {V<:Union{BigFloat,BigInterval}, T, U}
Performs the standard interval Newton method for a system of DynamicPolynomials.
Description
This function implements the same standard interval Newton method as above, with one key difference:
- The system is already in interval form (DynamicPolynomials with BigInterval coefficients)
- No input precision parameter is needed as the system is copied with the maximum precision of its coefficients
Arguments
sys::Vector
: A vector of DynamicPolynomials with BigInterval coefficientsinitial_point::Vector{V}
: Starting point (BigFloat or BigInterval)work_precision::Int64
: Precision for intermediate calculations. Default is 63target_precision::Int64
: Desired solution precision. Default is 2output_precision::Int64
: Precision for rounding final result. Default is 2max_nb_loops::Int64
: Maximum number of iterations. Default is 10refine::Int32
: Flag to control refinement behavior. Default is 1verbose::Cint
: Verbosity level. Default is 1
Returns
Same as above
interval_newton(f::DynamicPolynomials.Polynomial{T, U, BigInterval},
initial_point::V,
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10;
refine::Int32=Int32(1),
verbose::Cint=Cint(1)) where {V<:Union{BigFloat,BigInterval}, T, U}
Single equation version of interval Newton method for DynamicPolynomials.
Arguments
f
: A single DynamicPolynomial with BigInterval coefficientsinitial_point::V
: Starting point (BigFloat or BigInterval)- Other arguments as above
Returns
Same as above, but for a single equation
interval_newton(sys_ptr::Ptr{Cvoid},
initial_point::Vector{T},
work_precision::Int64=63,
target_precision::Int64=2,
output_precision::Int64=2,
max_nb_loops::Int64=10,
refine::Int32=Int32(1);
verbose::Cint=Cint(1)) where {T<:Union{BigFloat,BigInterval}}
Low-level version of the standard interval Newton method that works directly with C pointers.
Arguments
sys_ptr::Ptr{Cvoid}
: A C pointer to the polynomial system (typically obtained fromsys_get_C_ptr
)initial_point::Vector{T}
: Starting point (BigFloat or BigInterval)work_precision::Int64
: Precision for intermediate calculations. Default is 63target_precision::Int64
: Desired solution precision. Default is 2output_precision::Int64
: Precision for rounding final result. Default is 2max_nb_loops::Int64
: Maximum number of iterations. Default is 10refine::Int32
: Flag to control refinement behavior. Default is 1verbose::Cint
: Verbosity level. Default is 1
Returns
Tuple (status, solution, constants)
where:
status::Int
: Result code (0: singular Jacobian, 1: solution found, 2: max precision reached, 3: max iterations reached)solution::Vector{BigInterval}
: Computed solution intervalsconstants::Tuple
: Kantorovich test constants
Notes
- This is an internal function that should typically be called through the higher-level interfaces
- The system pointer must be properly initialized using
sys_get_C_ptr
Examples
# Using AbstractAlgebra
using AbstractAlgebra
R, (x,y) = polynomial_ring(QQ, ["x", "y"])
f = x^2 + y^2 - 1
g = x - y
sys = [f, g]
initial = [BigFloat(0.5), BigFloat(0.5)]
status, sol, consts = LACE.interval_newton(sys, initial)
# Using DynamicPolynomials
using DynamicPolynomials, MPFI
@polyvar x y
f = BigInterval(1)*x^2 + BigInterval(1)*y^2 - BigInterval(1.0, 1.1)
g = x - y
sys = [f, g]
initial = [BigInterval(0.5), BigInterval(0.5)]
status, sol, consts = LACE.interval_newton(sys, initial)
RationalUnivariateRepresentation.jl functions
RationalUnivariateRepresentation.zdim_parameterization
— Functionzdim_parameterization(system)
Computes a rational univariate representation of a system of polynomials.
Arguments
system
: an array of polynomials over the rationals. The input must form a zero-dimensional ideal, otherwise, an error is raised. Polynomials from these frontends are supported: AbstractAlgebra.jl, Nemo.jl.
Returns
The function zdim_parametrization
returns a single object rur
:
rur
: an array of polynomials, a rational univariate representation of the roots of the system; a polynomial is given by an array of its coefficients. Let $n$ be the number of indeterminates;rur
is an array $f(T), f_1(T), ..., f_n(T)$, such that ${(x_1,...,x_n) s.t. f(u) = 0, x_i = f_i(u) / f_0(u)}$, where $f_0$ is the derivative of the square-free part of $f$, is in bijection with the roots ofsystem
.
If the optional argument get_separating_element=true
is set, then the result is a tuple (rur
, sep
):
sep
: an array of coefficients $(a_1,..., a_n)$, which defines the linear form $T = a_1 x_1 + ... + a_n x_n$, that is the inverse of the bijection of the roots of thesystem
and the roots of $f(T)$.
Possible Options
The function zdim_parametrization
has the following optional arguments:
verbose::Bool
: a bool, whether to print progress statements or not. Default istrue
.get_separating_element::Bool
: a bool, whether to also return a separating element. Default isfalse
.composite::Int
: an Int, a power of two, the width of composite numbers in multi-modular computation. Default is4
.nn::Int32
: the bitsize of primes in multi-modular computation. Must be not exceed30
. Default is28
.parallelism::Symbol
: a symbol, parallelism mode. Available settings are:parallelism=:serial
: no multi-threading (default).parallelism=:multithreading
: use Julia threads. The number of Julia threads can be specified viathreads
option.
threads::Int
: an integer, the number of threads. Can only be used together withparallelism=:multithreading
. Default isBase.Threads.nthreads()
.search_strategy
: a Symbol, a search strategy. If the last variable is not separating, search for a separating form using the strategy. One of the following::current
(default): Original strategy.:random
: Use a1 x1 +...+ an xn, where ai is selected from [-B, B] for i = 1..n.:deterministic
: Try x1 + i x2 +...+ (i)^(n-1) xn for i = 1..n D (D-1) / 2, one by one.
Example
Using RUR.jl with AbstractAlgebra.jl:
julia> using AbstractAlgebra, RationalUnivariateRepresentation
julia> R, (x,y,z) = QQ["x","y","z"];
julia> system = [x*y*z - 1, x^2*y - 1, x*z - y];
julia> zdim_parameterization(system, verbose=false)
4-element Vector{Vector{Rational{BigInt}}}:
[-1, 0, 0, 0, 1]
[4]
[0, 4]
[0, 0, 0, 0, 4]
DiscriminantVariety.jl functions
DiscriminantVariety.discriminant_variety
— Functiondiscriminant_variety(sys, vars, params)
Computes a Discriminant Variety of system sys
from Q[params][vars]
.
RS.jl functions
PACE.isolate
— Functionisolate(args...; kwargs...)
A wrapper function for RS.rs_isolate
that isolates the real roots of a polynomial or polynomial system.
See also
PACE.RS.rs_isolate
Native functions
PACE.compute_system_solutions_from_rur
— Functioncompute_system_solutions_from_rur(sys::Vector,
target_precision::Int64=8;
verbose::Bool=false,
parallel::Bool=false,
refine::Cint=Cint(1))
Compute the solutions of a polynomial system using Rational Univariate Representation (RUR). The function first computes the RUR, then isolates the roots of the first polynomial, and finally finds the complete system solutions by evaluating the remaining RUR polynomials.
Arguments
sys::Vector
: A vector of polynomials representing the system of equationstarget_precision::Int64
: The target precision in bits for the mantissa of the solutions (default: 8)verbose::Bool
: Whether to print progress information during computation (default: false)parallel::Bool
: Whether to use parallel computation when possible (default: false)refine::Cint
: Controls the refinement strategy (default: 1)1
: Perform refinement using interval Newton method0
: No refinement
Returns
A vector of solutions to the polynomial system, where each solution is represented by intervals.
Notes
The computation follows these steps:
- Computes the RUR and separating element using
RationalUnivariateRepresentation.jl
- Isolates the roots of the first RUR polynomial using
RS.jl
- Evaluates the remaining RUR polynomials using
LACE.jl
to obtain the complete system solutions
compute_system_solutions_from_rur(rur::Vector{Vector{Rational{BigInt}}},
sep::Vector{Int64},
target_precision::Int64=8;
verbose::Bool=false,
parallel::Bool=false,
refine::Cint=Cint(1))
Compute the solutions of a polynomial system using a pre-computed Rational Univariate Representation (RUR). The function takes the RUR and separating element as input, computes the roots of the first RUR polynomial, and then finds the complete system solutions by evaluating the remaining RUR polynomials.
Arguments
rur::Vector{Vector{Rational{BigInt}}}
: The pre-computed RUR polynomialssep::Vector{Int64}
: The pre-computed separating elementtarget_precision::Int64
: The target precision in bits for the mantissa of the solutions (default: 8)verbose::Bool
: Whether to print progress information during computation (default: false)parallel::Bool
: Whether to use parallel computation when possible (default: false)refine::Cint
: Controls the refinement strategy (default: 1)1
: Perform refinement using interval Newton method0
: No refinement
Returns
A vector of solutions to the polynomial system, where each solution is represented by intervals.
Notes
This variant skips the RUR computation step and proceeds directly to:
- Isolating the roots of the first RUR polynomial using
RS.jl
- Evaluating the remaining RUR polynomials using
LACE.jl
to obtain the complete system solutions
PACE.partial_CAD
— Functionpartial_CAD(polynomials::Vector)
Compute a partial Cylindrical Algebraic Decomposition (CAD) for a system of polynomials and return sample rational points from each cell of the decomposition.
Arguments
polynomials::Vector
: A vector of polynomials defining the system
Returns
- A vector of sample points, where each point represents a cell in the CAD
Example
_, (x, y) = Nemo.polynomial_ring(Nemo.QQ, ["x", "y"])
# Define a system of polynomials
f1 = x - 3
f2 = x^2 + y^2 - 2
polys = [f1, f2]
# Get sample points from the CAD cells
sample_points = partial_CAD(polys)