Learn R Programming

R HIGHS Interface

Florian Schwendinger Updated: 2025-04-20

This repository contains an R interface to the HiGHS solver. The HiGHS solver, is a high-performance open-source solver for solving linear programming (LP), mixed-integer programming (MIP) and quadratic programming (QP) optimization problems.

1 Installation

The package can be installed from CRAN

install.packages("highs")

or GitLab.

remotes::install_gitlab("roigrp/solver/highs")

1.0.1 Using a preinstalled HiGHS library

It is possible to use a precompile HiGHS library by providing the system variable R_HIGHS_LIB_DIR. For example I used

mkdir build
cd build
cmake .. -DCMAKE_INSTALL_PREFIX=/Z/bin/highslib -DCMAKE_POSITION_INDEPENDENT_CODE:bool=ON -DSHARED:bool=OFF -DBUILD_TESTING:bool=OFF
make install

to install the HiGHS library to /Z/bin/highslib

Sys.setenv(R_HIGHS_LIB_DIR = "/Z/bin/highslib")
install.packages("highs")
# or 
# remotes::install_gitlab("roigrp/solver/highs")

2 Package

The highs package provides an API similar to Rglpk and a low level API to the HiGHS solver. For most users using highs_solve as shown below should be the best choice.

The package functions can be grouped into the following categories:

  1. The main function highs_solve.
  2. Object style wrappers for the model and the solver highs_model and highs_solver.
  3. Function highs_control to construct the control object for **highs_solve**and highs_solver.
  4. Low-level API wrapper functions to create and modify models hi_new_model and other functions starting with hi_model_.
  5. Low-level API wrapper functions to create and modify solvers hi_new_solver and other functions starting with hi_solver_.
  6. Functions to create example models example_model and solvers example_solver for the documentation examples.
  7. Function highs_available_solver_options to get the available solver options.
  8. Function highs_write_model to write the model to a file.
library("highs")

2.1 Examples

The the example models and solvers are included to have small examples available for the manual.

writeLines(ls("package:highs", pattern = "^example"))
#> example_model
#> example_solver

2.2 Low level model functions

The low-level model functions allow to create and modify models. More details and examples can be found in the manual.

writeLines(ls("package:highs", pattern = "^hi(|_new)_model"))
#> hi_model_get_ncons
#> hi_model_get_nvars
#> hi_model_set_constraint_matrix
#> hi_model_set_hessian
#> hi_model_set_lhs
#> hi_model_set_lower
#> hi_model_set_ncol
#> hi_model_set_nrow
#> hi_model_set_objective
#> hi_model_set_offset
#> hi_model_set_rhs
#> hi_model_set_sense
#> hi_model_set_upper
#> hi_model_set_vartype
#> hi_new_model

2.3 Low level solver functions

The low-level solver functions allow to create and modify solvers. More details and examples can be found in the manual.

writeLines(ls("package:highs", pattern = "^hi(|_new)_solver"))
#> hi_new_solver
#> hi_solver_add_cols
#> hi_solver_add_rows
#> hi_solver_add_vars
#> hi_solver_change_constraint_bounds
#> hi_solver_change_variable_bounds
#> hi_solver_clear
#> hi_solver_clear_model
#> hi_solver_clear_solver
#> hi_solver_get_bool_option
#> hi_solver_get_constraint_bounds
#> hi_solver_get_constraint_matrix
#> hi_solver_get_dbl_option
#> hi_solver_get_int_option
#> hi_solver_get_lp_costs
#> hi_solver_get_model
#> hi_solver_get_num_col
#> hi_solver_get_num_row
#> hi_solver_get_option
#> hi_solver_get_options
#> hi_solver_get_sense
#> hi_solver_get_str_option
#> hi_solver_get_variable_bounds
#> hi_solver_get_vartype
#> hi_solver_infinity
#> hi_solver_info
#> hi_solver_run
#> hi_solver_set_coeff
#> hi_solver_set_constraint_bounds
#> hi_solver_set_integrality
#> hi_solver_set_objective
#> hi_solver_set_offset
#> hi_solver_set_option
#> hi_solver_set_options
#> hi_solver_set_sense
#> hi_solver_set_variable_bounds
#> hi_solver_solution
#> hi_solver_status
#> hi_solver_status_message
#> hi_solver_write_basis
#> hi_solver_write_model

2.4 High level functions

The high level functions allow to work with models and solvers. More details and examples can be found in the manual.

args(highs_model)
#> function (Q = NULL, L, lower, upper, A = NULL, lhs = NULL, rhs = NULL, 
#>     types = rep.int(1L, length(L)), maximum = FALSE, offset = 0) 
#> NULL
args(highs_solver)
#> function (model, control = highs_control()) 
#> NULL
args(highs_control)
#> function (threads = 1L, time_limit = Inf, log_to_console = FALSE, 
#>     ...) 
#> NULL
args(highs_write_model)
#> function (model, file) 
#> NULL

2.5 Main function

The main function highs_solve.

library("highs")

args(highs_solve)
#> function (Q = NULL, L, lower, upper, A = NULL, lhs = NULL, rhs = NULL, 
#>     types = rep.int(1L, length(L)), maximum = FALSE, offset = 0, 
#>     control = highs_control()) 
#> NULL

2.6 LP

# Minimize
#  x_0 +  x_1 + 3
# Subject to
#                 x_1 <= 7
#  5 <=   x_0 + 2 x_1 <= 15
#  6 <= 3 x_0 + 2 x_1
#  0 <=   x_0         <= 4
#  1 <=           x_1
A <- rbind(c(0, 1), c(1, 2), c(3, 2))
s <- highs_solve(L = c(1.0, 1), lower = c(0, 1), upper = c(4, Inf),
                 A = A, lhs = c(-Inf, 5, 6), rhs = c(7, 15, Inf),
                 offset = 3)
str(s)
#> List of 6
#>  $ primal_solution: num [1:2] 0.5 2.25
#>  $ objective_value: num 5.75
#>  $ status         : int 7
#>  $ status_message : chr "Optimal"
#>  $ solver_msg     :List of 6
#>   ..$ value_valid: logi TRUE
#>   ..$ dual_valid : logi TRUE
#>   ..$ col_value  : num [1:2] 0.5 2.25
#>   ..$ col_dual   : num [1:2] 0 0
#>   ..$ row_value  : num [1:3] 2.25 5 6
#>   ..$ row_dual   : num [1:3] 0 0.25 0.25
#>  $ info           :List of 18
#>   ..$ valid                     : logi TRUE
#>   ..$ mip_node_count            : num -1
#>   ..$ simplex_iteration_count   : int 0
#>   ..$ ipm_iteration_count       : int 5
#>   ..$ qp_iteration_count        : int 0
#>   ..$ crossover_iteration_count : int 0
#>   ..$ primal_solution_status    : chr "Feasible"
#>   ..$ dual_solution_status      : chr "Feasible"
#>   ..$ basis_validity            : int 1
#>   ..$ objective_function_value  : num 5.75
#>   ..$ mip_dual_bound            : num 0
#>   ..$ mip_gap                   : num Inf
#>   ..$ num_primal_infeasibilities: int 0
#>   ..$ max_primal_infeasibility  : num 0
#>   ..$ sum_primal_infeasibilities: num 0
#>   ..$ num_dual_infeasibilities  : int 0
#>   ..$ max_dual_infeasibility    : num 0
#>   ..$ sum_dual_infeasibilities  : num 0

2.7 QP

# Minimize
#  0.5 x^2 - 2 x + y
# Subject to
#  x <= 3
zero <- .Machine$double.eps * 100
Q <- rbind(c(1, 0), c(0, zero))
L <- c(-2, 1)
A <- t(c(1, 0))

cntrl <- list(log_dev_level = 0L)
s <- highs_solve(Q = Q, L = L, A = A, lhs = 0, rhs = 3, control = cntrl)
str(s)
#> List of 6
#>  $ primal_solution: num [1:2] 3 0
#>  $ objective_value: num -6
#>  $ status         : int 10
#>  $ status_message : chr "Unbounded"
#>  $ solver_msg     :List of 6
#>   ..$ value_valid: logi TRUE
#>   ..$ dual_valid : logi TRUE
#>   ..$ col_value  : num [1:2] 3 0
#>   ..$ col_dual   : num [1:2] 0 1
#>   ..$ row_value  : num 3
#>   ..$ row_dual   : num -2
#>  $ info           :List of 18
#>   ..$ valid                     : logi TRUE
#>   ..$ mip_node_count            : num -1
#>   ..$ simplex_iteration_count   : int 1
#>   ..$ ipm_iteration_count       : int 0
#>   ..$ qp_iteration_count        : int 0
#>   ..$ crossover_iteration_count : int 0
#>   ..$ primal_solution_status    : chr "Feasible"
#>   ..$ dual_solution_status      : chr "Infeasible"
#>   ..$ basis_validity            : int 1
#>   ..$ objective_function_value  : num -6
#>   ..$ mip_dual_bound            : num 0
#>   ..$ mip_gap                   : num Inf
#>   ..$ num_primal_infeasibilities: int 0
#>   ..$ max_primal_infeasibility  : num 0
#>   ..$ sum_primal_infeasibilities: num 0
#>   ..$ num_dual_infeasibilities  : int 1
#>   ..$ max_dual_infeasibility    : num 1
#>   ..$ sum_dual_infeasibilities  : num 1

3 Additional information

3.1 Sparse matrices

The HiGHs C++ library internally supports the matrix formats csc (compressed sparse column matrix) and csr (compressed Sparse Row array). The highs package currently supports the following matrix classes:

  1. "matrix" dense matrices,
  2. "dgCMatrix" compressed sparse column matrix from the Matrix package,
  3. "dgRMatrix" compressed sparse row matrix from the Matrix package,
  4. "matrix.csc" compressed sparse column matrix from the SparseM package,
  5. "matrix.csr" compressed sparse row matrix from the SparseM package,
  6. "simple_triplet_matrix" coordinate format from the slam package.

If the constraint matrix A is provided as dgCMatrix, dgRMatrix, matrix.csc or matrix.csr the underlying data is directly passed to HiGHs otherwise it is first transformed into the csc format an afterwards passed to HiGHs

library("Matrix")

A <- rbind(c(0, 1), c(1, 2), c(3, 2))
csc <- as(A, "CsparseMatrix")  # dgCMatrix
s0 <- highs_solve(L = c(1.0, 1), lower = c(0, 1), upper = c(4, Inf),
                  A = csc, lhs = c(-Inf, 5, 6), rhs = c(7, 15, Inf),
                  offset = 3)

csr <- as(A, "RsparseMatrix")  # dgRMatrix
s1 <- highs_solve(L = c(1.0, 1), lower = c(0, 1), upper = c(4, Inf),
                  A = csr, lhs = c(-Inf, 5, 6), rhs = c(7, 15, Inf),
                  offset = 3)

library("SparseM")

csc <- as.matrix.csc(A)
s2 <- highs_solve(L = c(1.0, 1), lower = c(0, 1), upper = c(4, Inf),
                  A = csc, lhs = c(-Inf, 5, 6), rhs = c(7, 15, Inf),
                  offset = 3)

csr <- as.matrix.csr(A)
s3 <- highs_solve(L = c(1.0, 1), lower = c(0, 1), upper = c(4, Inf),
                  A = csr, lhs = c(-Inf, 5, 6), rhs = c(7, 15, Inf),
                  offset = 3)

library("slam")
stm <- as.simple_triplet_matrix(A)
s4 <- highs_solve(L = c(1.0, 1), lower = c(0, 1), upper = c(4, Inf),
                  A = stm, lhs = c(-Inf, 5, 6), rhs = c(7, 15, Inf),
                  offset = 3)

4 Options

The function highs_available_solver_options lists the available solver options

d <- highs_available_solver_options()
d[["option"]] <- sprintf("`%s`", d[["option"]])
knitr::kable(d, row.names = FALSE)
optiontype
presolvestring
solverstring
parallelstring
run_crossoverstring
time_limitdouble
read_solution_filestring
read_basis_filestring
write_model_filestring
solution_filestring
write_basis_filestring
random_seedinteger
rangingstring
infinite_costdouble
infinite_bounddouble
small_matrix_valuedouble
large_matrix_valuedouble
primal_feasibility_tolerancedouble
dual_feasibility_tolerancedouble
ipm_optimality_tolerancedouble
primal_residual_tolerancedouble
dual_residual_tolerancedouble
objective_bounddouble
objective_targetdouble
threadsinteger
user_bound_scaleinteger
user_cost_scaleinteger
highs_debug_levelinteger
highs_analysis_levelinteger
simplex_strategyinteger
simplex_scale_strategyinteger
simplex_crash_strategyinteger
simplex_dual_edge_weight_strategyinteger
simplex_primal_edge_weight_strategyinteger
simplex_iteration_limitinteger
simplex_update_limitinteger
simplex_min_concurrencyinteger
simplex_max_concurrencyinteger
log_filestring
write_model_to_filebool
write_presolved_model_to_filebool
write_solution_to_filebool
write_solution_styleinteger
glpsol_cost_row_locationinteger
write_presolved_model_filestring
output_flagbool
log_to_consolebool
timeless_logbool
ipm_iteration_limitinteger
pdlp_native_terminationbool
pdlp_scalingbool
pdlp_iteration_limitinteger
pdlp_e_restart_methodinteger
pdlp_d_gap_toldouble
qp_iteration_limitinteger
qp_nullspace_limitinteger
iis_strategyinteger
blend_multi_objectivesbool
log_dev_levelinteger
log_githashbool
solve_relaxationbool
allow_unbounded_or_infeasiblebool
use_implied_bounds_from_presolvebool
lp_presolve_requires_basis_postsolvebool
mps_parser_type_freebool
use_warm_startbool
keep_n_rowsinteger
cost_scale_factorinteger
allowed_matrix_scale_factorinteger
allowed_cost_scale_factorinteger
ipx_dualize_strategyinteger
simplex_dualize_strategyinteger
simplex_permute_strategyinteger
max_dual_simplex_cleanup_levelinteger
max_dual_simplex_phase1_cleanup_levelinteger
simplex_price_strategyinteger
simplex_unscaled_solution_strategyinteger
presolve_reduction_limitinteger
restart_presolve_reduction_limitinteger
presolve_substitution_maxfillininteger
presolve_rule_offinteger
presolve_rule_loggingbool
presolve_remove_slacksbool
simplex_initial_condition_checkbool
no_unnecessary_rebuild_refactorbool
simplex_initial_condition_tolerancedouble
rebuild_refactor_solution_error_tolerancedouble
dual_steepest_edge_weight_error_tolerancedouble
dual_steepest_edge_weight_log_error_thresholddouble
dual_simplex_cost_perturbation_multiplierdouble
primal_simplex_bound_perturbation_multiplierdouble
dual_simplex_pivot_growth_tolerancedouble
presolve_pivot_thresholddouble
factor_pivot_thresholddouble
factor_pivot_tolerancedouble
start_crossover_tolerancedouble
less_infeasible_DSE_checkbool
less_infeasible_DSE_choose_rowbool
use_original_HFactor_logicbool
run_centringbool
max_centring_stepsinteger
centring_ratio_tolerancedouble
icrashbool
icrash_dualizebool
icrash_strategystring
icrash_starting_weightdouble
icrash_iterationsinteger
icrash_approx_iterinteger
icrash_exactbool
icrash_breakpointsbool
mip_detect_symmetrybool
mip_allow_restartbool
mip_max_nodesinteger
mip_max_stall_nodesinteger
mip_max_start_nodesinteger
mip_max_leavesinteger
mip_max_improving_solsinteger
mip_lp_age_limitinteger
mip_pool_age_limitinteger
mip_pool_soft_limitinteger
mip_pscost_minreliableinteger
mip_min_cliquetable_entries_for_parallelisminteger
mip_report_levelinteger
mip_feasibility_tolerancedouble
mip_rel_gapdouble
mip_abs_gapdouble
mip_heuristic_effortdouble
mip_min_logging_intervaldouble
mip_heuristic_run_rinsbool
mip_heuristic_run_rensbool
mip_heuristic_run_root_reduced_costbool
mip_heuristic_run_zi_roundbool
mip_heuristic_run_shiftingbool
mip_improving_solution_savebool
mip_improving_solution_report_sparsebool
mip_improving_solution_filestring
mip_root_presolve_onlybool
mip_lifting_for_probinginteger

for additional information see the HiGHS homepage.

5 Status codes

HiGHS currently has the following status codes defined in HConst.h".

enumeratorstatusmessage
kNotset0"Not Set"
kLoadError1"Load error"
kModelError2"Model error"
kPresolveError3"Presolve error"
kSolveError4"Solve error"
kPostsolveError5"Postsolve error"
kModelEmpty6"Empty"
kOptimal7"Optimal"
kInfeasible8"Infeasible"
kUnboundedOrInfeasible9"Primal infeasible or unbounded"
kUnbounded10"Unbounded"
kObjectiveBound11"Bound on objective reached"
kObjectiveTarget12"Target for objective reached"
kTimeLimit13"Time limit reached"
kIterationLimit14"Iteration limit reached"
kUnknown15"Unknown"
kMin0"Not Set"
kMax15"Unknown"

Copy Link

Version

Install

install.packages('highs')

Monthly Downloads

2,601

Version

1.10.0-1

License

GPL (>= 2)

Maintainer

Florian Schwendinger

Last Published

April 24th, 2025

Functions in highs (1.10.0-1)

hi_reset_global_scheduler

Reset Global Scheduler
hi_model_set_offset

Set Offset for Highs Model
hi_model_set_objective

Set Objective for Highs Model
hi_new_model

Create new Highs Model
hi_solver_get_constraint_bounds

Get Constraint Bounds
hi_solver_set_coeff

Set a coefficient in the constraint matrix.
hi_solver_get_str_option

Get String Option Value
hi_solver_get_bool_option

Get Boolean Option Value
hi_new_solver

Create a new solver instance.
hi_solver_get_variable_bounds

Get Variable Bounds
hi_solver_get_vartype

Get Variable Types
hi_solver_add_vars

Add Variables to the Solver
hi_solver_add_rows

Add Constraints to Model
hi_solver_clear_solver

Clear Solver State
hi_solver_run

Run the Solver
hi_solver_get_solution

Get Solution
highs_available_solver_options

Available Solver Options
hi_solver_get_dbl_option

Get Double Option Value
hi_solver_get_constraint_matrix

Get Constraint Matrix
hi_solver_get_options

Get multiple HiGHS Solver Options
hi_solver_get_sense

Get the optimization sense of the solver instance.
hi_model_set_upper

Set Upper Bounds for a Highs Model
hi_solver_clear

Clear All Solver Data
hi_solver_clear_model

Clear Model Data
hi_solver_status_message

Get Solver Status Message
hi_solver_status

Get Solver Status
highs_write_model

Write a Highs Model to a File
hi_solver_set_constraint_bounds

Set constraint bounds for a given constraint.
highs_solver

Highs Solver
hi_solver_set_integrality

Set integrality for a variable in the solver.
highs_control

Highs Control
hi_solver_get_num_row

Get Number of Constraints
hi_solver_set_objective

Set the objective coefficient for a variable.
hi_solver_get_option

Get a HiGHS Solver Option
hi_solver_set_offset

Set the objective offset for the solver.
hi_solver_set_sense

Set the optimization sense of the solver instance.
hi_solver_set_variable_bounds

Set variable bounds for a given variable.
hi_model_set_vartype

Set Variable Types in a Highs Model
hi_solver_set_option

Set a HiGHS Solver Option
hi_solver_get_int_option

Get Integer Option Value
hi_solver_set_options

Set Multiple HiGHS Solver Options
highs_model

Create a Highs Model
hi_solver_change_constraint_bounds

Change Constraint Bounds
hi_solver_info

Get Solver Information
hi_solver_get_lp_costs

Get Objective Coefficients
hi_solver_get_num_col

Get Number of Variables
hi_solver_infinity

Get Solver Infinity Value
hi_solver_change_variable_bounds

Change Variable Bounds
highs_solve

Solve an Optimization Problems
hi_solver_write_basis

Write Basis to File
hi_solver_write_model

Write Model to File
hi_model_set_hessian

Set Hessian Matrix for Highs Model
hi_model_set_constraint_matrix

Set Constraint Matrix for Highs Model
hi_model_set_ncol

Sets the number of columns in the model
hi_model_set_lhs

Set Left Hand Side for a Highs Model
hi_model_set_nrow

Set the number of rows in the model
hi_model_get_ncons

Get Number of Constraints in a Model
hi_model_get_nvars

Get Number of Variables in a Highs Model
example_solver

Create a HiGHS Solver Object
example_model

Generate Example Optimization Models
hi_model_set_lower

Set Lower Bounds for Highs Model
hi_model_set_sense

Set the sense of the optimization model
hi_model_set_rhs

Set Right Hand Side for a Highs Model
hi_solver_add_cols

Add Variables to Model