Skip to content

test: add JET and fix issues#687

Open
oameye wants to merge 5 commits intoJuliaHomotopyContinuation:mainfrom
oameye:jet
Open

test: add JET and fix issues#687
oameye wants to merge 5 commits intoJuliaHomotopyContinuation:mainfrom
oameye:jet

Conversation

@oameye
Copy link
Copy Markdown
Contributor

@oameye oameye commented Mar 19, 2026

Fix 148 JET.jl static analysis reports (151 → 3)

Hey, I am looking to improve the TTFX of the package a bit. I thought a good start would be to add JET.jl to the test suite.
This PR used Claude during development to solve the JET reports. The overview below was also made with Claude.
It found a couple of serious bugs like:

function Base.mod(x::DoubleF64, y::DoubleF64)
    n = round(a / b)
    return (a - b * n)
end

However, since the current test suite didn't catch it. These methods can likely be removed?

Anyway let me know what you think.

Summary

Ran JET.report_package(HomotopyContinuation; target_modules=(HomotopyContinuation,)) and systematically fixed 148 of 151 reports. The remaining 3 are filtered as known (cos/sin/sqrt for IComplexF64 — implementing them causes +4s TTFX regression via @generated execute_instructions! recompilation).

Zero performance regression on precompilation, load time, TTFX, and steady-state runtime. Full test suite passes (2549/2549).

Highlights

Real bugs found (12 fixes)

JET uncovered genuine bugs that would crash at runtime — meaning these code paths are never exercised by the test suite and may be candidates for removal:

  • DoubleDouble.jl:351mod(x, y) used undefined a, b (copy-paste from divrem above)
  • interval_arithmetic.jl:50hash(::Interval, h) used undefined u instead of h
  • interval_arithmetic.jl:249inv(pow(x, -n)) where pow is never defined
  • symbolic.jl:789rethrow(e) but the catch variable is err
  • symbolic.jl:1102,1442map(var_groups) before var_groups is assigned (should be map(variable_groups))
  • intermediate_representation.jl:384throw(ExprError(...)) but ExprError is never defined
  • interpreted_homotopy.jl:221 — Acb evaluate! missing t parameter entirely
  • linear_algebra.jl:342 — unqualified SingularException (needs LA.SingularException)
  • utils.jl:484/498/513mul!(c::Float64, x::BigFloat) = mul!(x, c) — 2-arg form doesn't exist (infinite recursion)
  • norm.jl:20,27distance/norm fallbacks return a MethodError object instead of throwing it
  • DoubleDouble.jl:1032,1046DomainError() with no arguments

These all crash unconditionally when reached, proving the code paths have zero test coverage. In particular mod(::DoubleF64, ::DoubleF64), hash(::Interval), Interval^(-n), and the BigFloat mul!/add! 2-arg forms are completely broken — worth considering whether they should be removed rather than fixed.

Type narrowing (JET can't prove runtime invariants)

The bulk of fixes help JET verify types through patterns it can't reason about. For example, many mutable structs store Union{Nothing,T} fields that are lazily initialized:

# Before: JET reports no matching method `setprecision!(::Nothing, ...)`
if isnothing(H.eval_acb)
    H.eval_acb = interpreter(AcbRefVector, H.eval_ComplexF64)
end
setprecision!(H.eval_acb, prec)  # JET still sees Union{Nothing, Interpreter}

# After: type assertion tells JET (and the runtime) that Nothing is impossible here
if isnothing(H.eval_acb)
    H.eval_acb = interpreter(AcbRefVector, H.eval_ComplexF64)
end
I = H.eval_acb::Interpreter{AcbRefVector}  # narrows type; throws TypeError if wrong
setprecision!(I, prec)

For the same reason something() is used throughout — it unwraps Union{Nothing,T} and throws a clear ArgumentError("nothing") instead of a confusing MethodError:

# Before: JET reports no matching method `overlaps(::Nothing, ::Nothing)`
Arblib.overlaps(cert.I, certᵢ.I)

# After: something() unwraps or throws with a clear message
Arblib.overlaps(something(cert.I), something(certᵢ.I))

Other type narrowing patterns used:

  • Struct field types narrowed from Union{Nothing,Vector} to Vector where constructors always provide concrete values (subspace_homotopies.jl a_minus_b/offset)
  • Type assertions at union-split call sites (F::RandomizedSystem, cell::MixedCell, NewtonCache(...)::NewtonCache{...})
  • Type annotations on untyped function parameters throughout numerical_irreducible_decomposition.jl and witness_set.jl

Dispatch / missing method fixes

  • on_affine_chart(::AbstractSystem) had compile as positional instead of keyword — actual dispatch failure
  • 5 constructors (CompiledSystem, InterpretedSystem, CompiledHomotopy, InterpretedHomotopy, MixedHomotopy) didn't accept kwargs... passed by fixed()
  • Missing methods: variable_groups(::Homotopy), Homotopy(::InterpretedHomotopy), path_results(::AbstractVector{<:PathResult}), jacobian!(::InterpretedSystem, 5-arg), _variables(::ExpressionRef), horner(::Expression, ::Variable)
  • MatrixWorkspace constructor refactored for type-stability (extracted helper to avoid Union{Matrix, StructArray} in constructor body)

Filtered reports (3)

cos/sin/sqrt for IComplexF64 are not implemented. Adding them (whether on Base.cos or on the package-internal op_cos) triggers recompilation of the @generated execute_instructions! interpreter for all IComplex tape types, causing +4s TTFX on first certify(). These only matter when certifying systems with transcendental functions. The test filters them with a comment explaining the trade-off.

Test integration

Added test/jet_test.jl as part of the regular test suite (JET added as test dependency). Filters the 3 known IComplexF64 reports and asserts 0 remaining.

Benchmarks (vs main)

Metric Verdict
Precompilation No change (9.49s vs 9.51s)
Load time No change (1.70s vs 1.71s)
TTFX solve No change (32.7s vs 33.0s)
TTFX certify No change (20.6s vs 20.3s)
Runtime (cyclic-5, certify, param) No change

Benchmarked sequentially (no concurrency), alternating branches, 3 runs each.

  • Precompilation: @elapsed Pkg.precompile() after clearing ~/.julia/compiled/*HomotopyContinuation*

  • Load time: @elapsed using HomotopyContinuation with cache present

  • TTFX: @time on first call (includes JIT):

    using HomotopyContinuation
    @var x y
    F = System([x^2 + y^2 - 1, x - y])
    @time solve(F; show_progress=false)                    # TTFX solve
    F2 = System([x^2 + y^2 - 1, x*y - 1])
    res = solve(F2; show_progress=false)
    @time certify(F2, res; show_progress=false)            # TTFX certify
  • Runtime: @benchmark after full warmup:

    using HomotopyContinuation, BenchmarkTools
    @var x y z[1:5] a b
    
    # cyclic-5
    n = 5
    polys = [sum(prod(z[mod1(j+k,n)] for j in 1:i) for k in 1:n) - (-1)^(i-1)*i for i in 1:n-1]
    push!(polys, prod(z) - 1); F_c = System(polys)
    solve(F_c; show_progress=false)  # warmup
    @benchmark solve($F_c; show_progress=false)
    
    # certify
    F2 = System([x^2 + y^2 - 1, x*y - 1])
    r = solve(F2; show_progress=false); certify(F2, r; show_progress=false)
    @benchmark certify($F2, $r; show_progress=false)
    
    # parameter homotopy
    G = System([x^2 - a, y^2 - b]; parameters=[a, b])
    solve(G; target_parameters=[1.0, 1.0], show_progress=false)
    @benchmark solve($G; target_parameters=[2.0, 3.0], show_progress=false)
    Benchmark main (median) JET (median) main (allocs) JET (allocs)
    cyclic-5 19.59 ms 19.73 ms 23732 23732
    certify 362.9 μs 373.4 μs 3127 3142
    param homotopy 547.0 μs 546.9 μs 4776 4775

    No runtime regression. Identical allocation counts (±15 allocs on certify from the split Arblib method).

julia> versioninfo()
Julia Version 1.12.5
Commit 5fe89b8ddc1 (2026-02-09 16:05 UTC)
Build Info:
  Official https://julialang.org release
Platform Info:
  OS: Linux (x86_64-linux-gnu)
  CPU: 12 × AMD Ryzen 5 5600X 6-Core Processor
  WORD_SIZE: 64
  LLVM: libLLVM-18.1.7 (ORCJIT, znver3)
  GC: Built with stock GC
Threads: 1 default, 1 interactive, 1 GC (on 12 virtual cores)

@PBrdng
Copy link
Copy Markdown
Collaborator

PBrdng commented Mar 19, 2026

Thanks! I pushed after running JuliaFormatter (we still use an old version, should update it at some point...).

@oameye
Copy link
Copy Markdown
Contributor Author

oameye commented Mar 19, 2026

btw I think it would also be good to add AQUA.jl. However this will likely require bigger changes.

@oameye
Copy link
Copy Markdown
Contributor Author

oameye commented Mar 19, 2026

Seems that the JET test failed on lts. I didn't test on 1.10, and since the compiler changes between Julia version, it finds different issues. Will have a look later. Feel free to review already the other changes.

@oameye
Copy link
Copy Markdown
Contributor Author

oameye commented Mar 19, 2026

On LTS it report a false positive for the macro @var. Lts also only support an older version of JET. So I made the JET test only run on v1.12 or higher. Can you run the formatter again :)

@saschatimme
Copy link
Copy Markdown
Member

saschatimme commented Mar 19, 2026

Since tests take absolutely ages, probably best to split this into a separate check

@PBrdng
Copy link
Copy Markdown
Collaborator

PBrdng commented Mar 19, 2026

Since tests take absolutely ages, probably best to split this into a separate check

The test took very long already before this PR. It doesnt seem to spend significantly more time. But yes, splitting it into a separate check makes sense to keep the overview.

@oameye
Copy link
Copy Markdown
Contributor Author

oameye commented Mar 19, 2026

Yeah good idea to parallelize :)

@PBrdng
Copy link
Copy Markdown
Collaborator

PBrdng commented Mar 20, 2026

Yeah good idea to parallelize :)

I can take care of that in April (I dont have time before :()

@oameye
Copy link
Copy Markdown
Contributor Author

oameye commented Mar 21, 2026

I can take care of that in April (I dont have time before :()

Now worries. I am on a PhD thesis deadline until the end of the month. Will try to finish this afterwards :)

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

3 participants