API Reference

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

Main Functions

HPRQP.build_from_mpsFunction
build_from_mps(filename::String; verbose::Bool=true)

Build a QP model from an MPS file.

This function reads a QP problem from an MPS file and returns a CPU-based model that can be solved with optimize() or solve().

Arguments

  • filename::String: Path to the MPS file
  • verbose::Bool: Whether to print progress information (default: true)

Returns

  • QP_info_cpu: QP model ready to be solved

Example

using HPRQP

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

See also: build_from_QAbc, build_from_mat, optimize

source
HPRQP.build_from_QAbcFunction
build_from_QAbc(Q, c, A, AL, AU, l, u, obj_constant=0.0)

Build a QP model from matrix form.

This function creates a QP problem from the standard form: min 0.5 <x,Qx> + <c,x> + obj_constant s.t. AL <= Ax <= AU l <= x <= u

Accepts both sparse and dense matrices for Q and A. Dense matrices will be automatically converted to sparse format for efficient computation.

Arguments

  • Q::Union{SparseMatrixCSC, Matrix{Float64}}: Quadratic objective matrix (n × n). Can be sparse or dense.
  • c::Vector{Float64}: Linear objective coefficients (length n)
  • A::Union{SparseMatrixCSC, Matrix{Float64}}: Constraint matrix (m × n). Can be sparse or dense.
  • 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

  • QP_info_cpu: QP model ready to be solved

Example

using SparseArrays, HPRQP

# Example 1: Sparse matrices
Q = sparse([2.0 0.0; 0.0 2.0])
c = [-3.0, -5.0]
A = sparse([-1.0 -2.0; -3.0 -1.0])
AL = [-10.0, -12.0]
AU = [Inf, Inf]
l = [0.0, 0.0]
u = [Inf, Inf]

model = build_from_QAbc(Q, c, A, AL, AU, l, u)
params = HPRQP_parameters()
result = optimize(model, params)

# Example 2: Dense matrices (automatically converted)
n = 10
Q = zeros(n, n)  # Empty or dense Q matrix
Q[1,1] = 2.0
c = ones(n)
A = ones(5, n)  # Dense constraint matrix
AL = -Inf * ones(5)
AU = ones(5)
l = zeros(n)
u = ones(n)
model = build_from_QAbc(Q, c, A, AL, AU, l, u)

See also: build_from_mps, build_from_mat, optimize

source
HPRQP.build_from_matFunction
build_from_mat(filename::String; problem_type::String="QAP", lambda::Float64=1.0, verbose::Bool=true)

Build a QP model from a MAT file (for QAP or LASSO problems).

This function reads a QAP (Quadratic Assignment Problem) or LASSO problem from a .mat file. Note: This function stores metadata that will be used to create operator-based models during solve().

Arguments

  • filename::String: Path to the MAT file
  • problem_type::String: Type of problem - "QAP" or "LASSO" (default: "QAP")
  • lambda::Float64: Regularization parameter for LASSO (default: 1.0)
  • verbose::Bool: Whether to print progress information (default: true)

Returns

  • Tuple: (metadatadict, problemtype) containing problem data

Example

using HPRQP

model_info, prob_type = build_from_mat("qap_problem.mat", problem_type="QAP")
params = HPRQP_parameters()
params.problem_type = prob_type
result = optimize((model_info, prob_type), params)

See also: build_from_mps, build_from_QAbc, optimize

source
HPRQP.build_from_ABSTFunction
build_from_ABST(A, B, S, T; verbose::Bool=true)

Build a QP model for Quadratic Assignment Problem (QAP) from matrices A, B, S, T.

This function creates a QAP problem in the standard form used by HPR-QP: min <vec(X), Qvec(X)> s.t. Xe = e, X'*e = e (doubly stochastic constraints) X >= 0

Where Q(X) = 2(AXB - SX - X*T) is represented as a matrix-free operator using the CUSTOMQOPERATOR API.

Arguments

  • A::Matrix{Float64}: Distance matrix for facility locations (n × n)
  • B::Matrix{Float64}: Flow matrix between facilities (n × n)
  • S::Matrix{Float64}: Linear term for rows (n × n)
  • T::Matrix{Float64}: Linear term for columns (n × n)
  • verbose::Bool: Whether to print progress information (default: true)

Returns

  • QP_info_cpu: QP model ready to be solved with operator-based Q

Example

using HPRQP

# Define QAP data matrices
n = 10
A = rand(n, n)
B = rand(n, n)
S = zeros(n, n)
T = zeros(n, n)

model = build_from_ABST(A, B, S, T)
params = HPRQP_parameters()
result = optimize(model, params)

See also: build_from_Ab_lambda, build_from_mat, build_from_QAbc, optimize

source
HPRQP.build_from_Ab_lambdaFunction
build_from_Ab_lambda(A, b, lambda; verbose::Bool=true)

Build a QP model for LASSO regression from data matrix A, target vector b, and regularization λ.

This function creates a LASSO problem in the standard form: min 0.5 ||A*x - b||₂² + λ ||x||₁

Which is reformulated as a QP with operator-based Q: min 0.5 <x, Q*x> + <c, x> + constant s.t. (no constraints on x, handled via proximal operator for L1)

Where Q = A'*A is represented as a matrix-free operator using the CUSTOMQOPERATOR API.

Arguments

  • A::SparseMatrixCSC{Float64}: Data matrix (m × n)
  • b::Vector{Float64}: Target vector (length m)
  • lambda::Float64: Regularization parameter (must be positive)
  • verbose::Bool: Whether to print progress information (default: true)

Returns

  • QP_info_cpu: QP model ready to be solved with operator-based Q

Example

using HPRQP, SparseArrays

# Define LASSO data
m, n = 100, 50
A = sprandn(m, n, 0.1)
b = randn(m)
lambda = 0.01 * norm(A' * b, Inf)

model = build_from_Ab_lambda(A, b, lambda)
params = HPRQP_parameters()
result = optimize(model, params)

See also: build_from_ABST, build_from_mat, build_from_QAbc, optimize

source
HPRQP.optimizeFunction
optimize(model::QP_info_cpu, params::HPRQP_parameters)

Solve a QP model with optional warm-up phase.

This is the main entry point for solving QP problems. It handles:

  1. Optional warm-up phase to avoid JIT compilation overhead
  2. Calls solve() which does scaling, GPU transfer, and optimization

Arguments

  • model::QP_info_cpu: QP model built from buildfrommps(), buildfromQAbc(), etc.
  • params::HPRQP_parameters: Solver parameters

Returns

  • HPRQP_results: Solution results

Example

using HPRQP

model = build_from_mps("problem.mps")
params = HPRQP_parameters()
params.stoptol = 1e-8
params.warm_up = true
result = optimize(model, params)

See also: build_from_mps, build_from_QAbc

source
HPRQP.OptimizerType
Optimizer()

Create a new HPRQP Optimizer object.

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

Example

using JuMP, HPRQP
model = JuMP.Model(HPRQP.Optimizer)
set_optimizer_attribute(model, "stoptol", 1e-6)
set_optimizer_attribute(model, "use_gpu", true)
set_optimizer_attribute(model, "device_number", 0)
source

Main Types

HPRQP_parameters

Solver parameters struct. Create with default values using HPRQP_parameters().

Key Parameters:

  • stoptol::Float64: Convergence tolerance (default: 1e-6)
  • max_iter::Int: Maximum iterations (default: 1000000)
  • use_gpu::Bool: Enable GPU acceleration (default: true)
  • verbose::Bool: Print iteration logs (default: true)
  • warm_up::Bool: Run warmup solve to eliminate JIT overhead (default: false)
  • time_limit::Float64: Maximum solve time in seconds (default: Inf)

See Parameters Guide for complete list.

HPRQP_results

Solution results struct returned by optimize().

Main Fields:

  • status::String: Solution status ("OPTIMAL", "MAXITER", "TIMELIMIT")
  • x::Vector{Float64}: Primal solution vector
  • y::Vector{Float64}: Dual solution for constraints
  • z::Vector{Float64}: Dual solution for bounds
  • primal_obj::Float64: Primal objective value
  • dual_obj::Float64: Dual objective value
  • iter::Int: Total iterations
  • time::Float64: Solve time in seconds

See Output & Results Guide for complete list.

Q Operators

The following Q operator types are available for specialized problems:

  • LASSO_Q_operator_cpu: For LASSO (L1-regularized least squares) problems
  • QAP_Q_operator_cpu: For quadratic assignment problems
  • Sparse matrices (SparseMatrixCSC): For general quadratic programming

See the Q Operators Guide for detailed information.

Quick Reference

Solving Problems

Direct API (Matrix Form):

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

# Step 2: Set parameters
params = HPRQP_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 = HPRQP_parameters()

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

JuMP:

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

Common Parameter Settings

params = HPRQP_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