Miscellaneous

Storing solutions and parameters

HomotopyContinuation.write_solutionsFunction
write_solutions(filename, solutions)

Stores the as a plain text file onto disk. The storage format is as follows. The first line indicates the number of solutions stored followed by a blank line. Then the solutions are stored where each solution is separated by a blank line. Note that the solutions are always considered as complex numbers. See read_solutions for reading the solution file back.

Note that this is the same file format as used in Bertini.

Example

julia> write_solutions("solutions.txt", [[1, 1], [-1, 2]])

shell> cat solutions.txt
2

1 0
1 0

-1 0
2 0

julia> read_solutions("solutions.txt")
2-element Array{Array{Complex{Int64},1},1}:
 [1 + 0im, 1 + 0im]
 [-1 + 0im, 2 + 0im]
source
HomotopyContinuation.read_solutionsFunction
read_solutions(filename)

Read the solutions stored with write_solutions.

Example

julia> write_solutions("solutions.txt", [[1, 1], [-1, 2]])

julia> read_solutions("solutions.txt")
2-element Array{Array{Complex{Int64},1},1}:
 [1 + 0im, 1 + 0im]
 [-1 + 0im, 2 + 0im]
source
HomotopyContinuation.write_parametersFunction
write_parameters(filename, parameters)

Stores the parameters as a plain text file onto disk. The storage format is as follows. The first line indicates the number of parameter values stored followed by a blank line. Then the parameter values are stored.

See read_parameters

Note that this is the same file format as used in Bertini.

Example

julia> write_parameters("parameters.txt", [2.0, -3.2 + 2im])

shell> cat parameters.txt
2

2.0 0.0
-3.2 2.0

julia> read_parameters("parameters.txt")
2-element Array{Complex{Float64},1}:
  2.0 + 0.0im
 -3.2 + 2.0im
source
HomotopyContinuation.read_parametersFunction
read_parameters(filename)

Read the parameters stored with write_parameters.

Example

julia> write_parameters("parameters.txt", [2.0, -3.2 + 2im])

julia> read_parameters("parameters.txt")
2-element Array{Complex{Float64},1}:
  2.0 + 0.0im
 -3.2 + 2.0im
source

Newton's method

HomotopyContinuation.newtonFunction
newton(
    F::AbstractSystem,
    x₀::AbstractVector,
    p = nothing,
    norm::AbstractNorm = InfNorm(),
    cache::NewtonCache = NewtonCache(F);
    options...
)

An implemenetation of a local Newton's method with various options to specify convergence criteria. Returns a NewtonResult. The computations are always performed in complex arithmetic with double precision, i.e., using Complex{Float64}. The optional cache argument pre-allocates the necessary memory. This is useful if the method is called repeatedly.

Options

  • atol::Float64 = 1e-8: The method is declared converged if norm(xᵢ₊₁ - xᵢ) < atol.
  • rtol::Float64 = atol: The method is declared converged if norm(xᵢ₊₁ - xᵢ) < max(atol, rtol * norm(x₀)).
  • max_iters::Int = 20: The maximal number of iterations.
  • extended_precision::Bool = false: An optional use of extended precision for the evaluation of F(x). This can increase the achievable accuracy.
  • contraction_factor::Float64 = 1.0: The Newton updates have to satisfy $|xᵢ₊₁ - xᵢ| < a^{2^{i-1}}|x₁ - x₀|$ for $i ≥ 1$ where $a$ is contraction_factor.
  • min_contraction_iters::Int = typemax(Int): The minimal number of iterations the contraction_factor has to be satisfied. If after min_contraction_iters many iterations the contraction factor is not satisfied the step is accepted anyway.
  • max_abs_norm_first_update::Float64 = Inf: The initial guess x₀ is rejected if norm(x₁ - x₀) > max_abs_norm_first_update
  • max_rel_norm_first_update::Float64 = max_abs_norm_first_update: The initial guess x₀ is rejected if norm(x₁ - x₀) > max_rel_norm_first_update * norm(x₀)
source
HomotopyContinuation.NewtonResultType
NewtonResult

Result returned by newton.

Fields

  • return_code::Symbol: Can be :success, :rejected or :max_iters.
  • x::Vector{ComplexF64}: The last value obtained.
  • accuracy::Float64: Estimate of the absolute distance of x to a true zero.
  • iters::Int Number of iterations performed.
  • contraction_ratio::Float64: The value |xᵢ - xᵢ₋₁| / |xᵢ₋₁ - xᵢ₋₂|.
source

Norms

HomotopyContinuation.AbstractNormType
AbstractNorm

An AbstractNorm represents any norm of a vector space. All norms are callable. norm(x) computes the norm of x and norm(x,y) computes the distance norm(x - y).

source
HomotopyContinuation.WeightedNormType
WeightedNorm(d::AbstractVector, norm::AbstractNorm; options...)
WeightedNorm(norm::AbstractNorm, n::Integer; options...)
WeightedNorm(norm::AbstractNorm, x::AbstractVector; options...)

A WeightedNorm represents a weighted variant of norm norm of a n-dimensional vector space.A norm||x||is weighted by introducing a vector of additional weightsdsuch that the new norm is||D⁻¹x||whereDis the diagonal matrix with diagonald. The WeightedNorm is desigened to change the weights dynamically by using init!(::WeightedNorm, x) and update!(::WeightedNorm, x). The weights are there constructed such that $||D⁻¹x|| ≈ 1.0$. The weights can be accessed and changed by indexing.

Options

  • scale_min = sqrt(eps()): The minimal size of dᵢ is scale_min time the (weighted) norm of x.
  • scale_abs_min = min(scale_min^2, 200 * sqrt(eps())): The absolute minimal size of dᵢ.
  • scale_max = 1.0 / eps() / sqrt(2): The absolute maximal size of dᵢ
source
LinearAlgebra.normMethod
norm(u, norm::AbstractNorm)

Compute the norm ||u|| with respect to the given norm norm.

source
HomotopyContinuation.update!Method
update!(w::WeightedNorm, x)

Update the weighted norm w for x, this will interpolate between the previous weights and the norm of x.

source

Unique points, group actions and multiplicities

HomotopyContinuation.UniquePointsType
UniquePoints{T, Id, M}

A data structure for assessing quickly whether a point is close to an indexed point as determined by the given distances function M. The distance function has to be a metric. The indexed points are only stored by their identifiers Id.

UniquePoints(v::AbstractVector{T}, id::Id;
                metric = EuclideanNorm(),
                group_actions = nothing)

Initialize the data structure. This does not initialize the data structure with the point.

Example

x = randn(ComplexF64, 4)
permutation(x) = ([x[2]; x[1]; x[3]; x[4]],)
group_actions = GroupActions(permutation)
X = group_actions(x)

# without group actions
unique_points = UniquePoints(x, 1)
HC.add!.(unique_points, X, 1:length(X), 1e-5)
length(unique_points) # 2

unique_points = UniquePoints(x, 1, group_actions = group_actions)
HC.add!.(unique_points, X, 1:length(X), 1e-5)
length(unique_points) # 1
source
HomotopyContinuation.search_in_radiusMethod
search_in_radius(unique_points, v, tol)

Search whether unique_points contains a point p with distances at most tol from v. Returns nothing if no point exists, otherwise the identifier of p is returned.

source
HomotopyContinuation.add!Method
add!(unique_points, v, id; atol = 1e-14, rtol = sqrt(eps()))
add!(unique_points, v, id, atol)

Search whether unique_points contains a point p with distances at most max(atol, norm(v)rtol) from v. If this is the case the identifier of p and false is returned. Otherwise (id, true) is returned.

source
HomotopyContinuation.multiplicitiesFunction
multiplicities(vectors; metric = EuclideanNorm(), atol = 1e-14, rtol = 1e-8, kwargs...)

Returns a Vector{Vector{Int}} v. Each vector w in 'v' contains all indices i,j such that w[i] and w[j] have distance at most max(atol, rtol * metric(0,w[i])). The remaining kwargs are things that can be passed to UniquePoints.

julia> multiplicities([[1,0.5], [1,0.5], [1,1]])
[[1,2]]

This is the same as

multiplicities([[1,0.5], [1,0.5], [1,1]]; distance=(x,y) -> LinearAlgebra.norm(x-y))

Here is an example for using group actions.

julia> X = [[1, 2, 3, 4], [2,1,3,4], [1,2,4,3], [2,1,4,3]]
julia> permutation(x) = [x[2], x[1], x[3], x[4]]
julia> m = multiplicities(X, group_action = permutation)
[[1,2], [3,4]]
source
HomotopyContinuation.unique_pointsFunction
unique_points(vectors; metric = EuclideanNorm(), atol = 1e-14, rtol = 1e-8, kwargs...)

Returns all elements in vector for which two elements have distance at most max(atol, rtol * metric(0,w[i])). Note that the output can depend on the order of elements in vectors. The remaining kwargs are things that can be passed to UniquePoints.

source

Debugging

HomotopyContinuation.path_infoFunction
path_info(tracker::Tracker, x₀, t₁ = 1.0, t₀ = 0.0; debug::Bool = false, kwargs...)

Track a path using the given tracker and start value x₀. This returns a struct containing detailed information about the tracked path.

source