API Reference

Complete API documentation for HPRLP. For detailed guides, see:

Main Functions

HPRLP.build_from_mpsFunction
build_from_mps(filename; verbose=true)

Build an LP model from an MPS file.

Arguments

  • filename::String: Path to the .mps file
  • verbose::Bool: Enable verbose output (default: true)

Returns

  • LP_info_cpu: LP model ready to be solved

Example

using HPRLP

model = build_from_mps("problem.mps")
params = HPRLP_parameters()
result = optimize(model, params)

See also: build_from_Abc, optimize

source
HPRLP.build_from_AbcFunction
build_from_Abc(A, c, AL, AU, l, u, obj_constant=0.0)

Build an LP model from matrix form.

Arguments

  • A::Union{SparseMatrixCSC, Matrix}: Constraint matrix (m × n). Dense matrices will be automatically converted to sparse format with a warning.
  • c::Vector{Float64}: Objective coefficients (length n)
  • AL::Vector{Float64}: Lower bounds for constraints Ax (length m)
  • AU::Vector{Float64}: Upper bounds for constraints Ax (length m)
  • l::Vector{Float64}: Lower bounds for variables x (length n)
  • u::Vector{Float64}: Upper bounds for variables x (length n)
  • obj_constant::Float64: Constant term in objective function (default: 0.0)

Returns

  • LP_info_cpu: LP model ready to be solved

Example

using SparseArrays, HPRLP

A = sparse([1.0 2.0; 3.0 1.0])
c = [-3.0, -5.0]
AL = [-Inf, -Inf]
AU = [10.0, 12.0]
l = [0.0, 0.0]
u = [Inf, Inf]

model = build_from_Abc(A, c, AL, AU, l, u)
params = HPRLP_parameters()
result = optimize(model, params)

See also: build_from_mps, optimize

source
HPRLP.optimizeFunction
optimize(model::LP_info_cpu, params::HPRLP_parameters)

Optimize a linear program using the HPR-LP algorithm.

This function handles GPU transfer, scaling, and optional warmup internally based on the parameters.

Arguments

  • model::LP_info_cpu: LP model built from build_from_mps or build_from_Abc
  • params::HPRLP_parameters: Solver parameters

Returns

  • HPRLP_results: Solution results including objective value, solution vector, and convergence info

Example

using HPRLP

model = build_from_mps("problem.mps")
params = HPRLP_parameters()
params.stoptol = 1e-6
params.use_gpu = true
params.warm_up = true

result = optimize(model, params)
println("Status: ", result.status)
println("Objective: ", result.primal_obj)

See also: build_from_mps, build_from_Abc, HPRLP_parameters

source
HPRLP.OptimizerType
Optimizer()

Create a new HPRLP Optimizer object.

Set optimizer attributes using MOI.RawOptimizerAttribute or JuMP.set_optimizer_attribute.

Example

using JuMP, HPRLP
model = JuMP.Model(HPRLP.Optimizer)
set_optimizer_attribute(model, "stoptol", 1e-4)
set_optimizer_attribute(model, "use_gpu", true)
source

Types

HPRLP.HPRLP_parametersType
HPRLP_parameters

Parameters for the HPR-LP solver.

Fields

  • stoptol::Float64: Stopping tolerance (default: 1e-4)
  • max_iter::Int: Maximum number of iterations (default: typemax(Int32))
  • time_limit::Float64: Time limit in seconds (default: 3600.0)
  • check_iter::Int: Interval for residual checks (default: 150)
  • use_Ruiz_scaling::Bool: Enable Ruiz scaling (default: true)
  • use_Pock_Chambolle_scaling::Bool: Enable Pock-Chambolle scaling (default: true)
  • use_bc_scaling::Bool: Enable b/c scaling (default: true)
  • use_gpu::Bool: Use GPU acceleration (default: true)
  • device_number::Int: GPU device number (default: 0)
  • warm_up::Bool: Enable warm-up phase (default: true)
  • print_frequency::Int: Print log every N iterations, -1 for auto (default: -1)
  • verbose::Bool: Enable verbose output (default: true)
  • initial_x::Union{Vector{Float64},Nothing}: Initial primal solution (default: nothing)
  • initial_y::Union{Vector{Float64},Nothing}: Initial dual solution (default: nothing)
  • auto_save::Bool: Automatically save best x, y, and sigma during optimization (default: false)
  • save_filename::String: Filename for auto-save HDF5 file (default: "hprlp_autosave.h5")

Example

params = HPRLP_parameters()
params.stoptol = 1e-6
params.use_gpu = false
params.verbose = false
params.auto_save = true
params.save_filename = "my_optimization.h5"
# Provide initial point
params.initial_x = x0
params.initial_y = y0
source
HPRLP.HPRLP_resultsType
HPRLP_results

Results from the HPR-LP solver.

Fields

  • iter::Int: Total number of iterations
  • iter_4::Int: Iterations to reach 1e-4 accuracy
  • iter_6::Int: Iterations to reach 1e-6 accuracy
  • iter_8::Int: Iterations to reach 1e-8 accuracy
  • time::Float64: Total solve time in seconds
  • time_4::Float64: Time to reach 1e-4 accuracy
  • time_6::Float64: Time to reach 1e-6 accuracy
  • time_8::Float64: Time to reach 1e-8 accuracy
  • primal_obj::Float64: Primal objective value
  • dual_obj::Float64: Dual objective value
  • primal_residual::Float64: Final primal feasibility residual
  • dual_residual::Float64: Final dual feasibility residual
  • relative_duality_gap::Float64: Final relative duality gap
  • x::Vector{Float64}: Primal solution vector
  • status::String: Termination status ("Optimal", "TimeLimit", "IterationLimit")

Example

model = build_from_mps("problem.mps")
params = HPRLP_parameters()
results = solve(model, params)
println("Status: ", results.status)
println("Objective: ", results.primal_obj)
println("Time: ", results.time, " seconds")
source

Quick Reference

Solving Problems

Direct API (Matrix Form):

# Step 1: Build the model
model = build_from_Abc(A, c, AL, AU, l, u)

# Step 2: Set parameters
params = HPRLP_parameters()
params.stoptol = 1e-6

# Step 3: Optimize
result = optimize(model, params)

MPS Files:

# Step 1: Build the model from file
model = build_from_mps("problem.mps")

# Step 2: Set parameters
params = HPRLP_parameters()

# Step 3: Optimize
result = optimize(model, params)

JuMP:

model = Model(HPRLP.Optimizer)
# ... add variables and constraints ...
optimize!(model)

Common Parameter Settings

params = HPRLP_parameters()
params.stoptol = 1e-6           # Convergence tolerance
params.use_gpu = true           # Enable GPU
params.verbose = false          # Silent mode
params.time_limit = 3600        # Time limit (seconds)
params.warm_up = true           # Enable warmup for accurate timing
params.initial_x = x0           # Initial primal solution
params.initial_y = y0           # Initial dual solution
params.auto_save = true         # Auto-save best solution
params.save_filename = "opt.h5" # HDF5 file for auto-save

Accessing Results

result.status         # "OPTIMAL", "MAX_ITER", or "TIME_LIMIT"
result.primal_obj     # Primal objective value
result.x              # Primal solution vector
result.y              # Dual solution vector (constraints)
result.z              # Dual solution vector (bounds)
result.iter           # Total iterations
result.time           # Solve time (seconds)
result.residuals      # Final residual (max of primal, dual, gap)

See Also