Skip to content

Commit

Permalink
Merge pull request #7 from mckib2/embed-highs
Browse files Browse the repository at this point in the history
MPS writer and cleanup
  • Loading branch information
mdhaber committed Apr 18, 2020
2 parents 1a18376 + 25521e5 commit 2aa31fa
Show file tree
Hide file tree
Showing 6 changed files with 106 additions and 41 deletions.
6 changes: 3 additions & 3 deletions scipy/optimize/_highs/PYREADME.rst
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@ Cython wrappers over `HiGHS <https://github.com/ERGO-Code/HiGHS>`_. This projec

The package is available on pypi as `scikit-highs <https://pypi.org/project/scikit-highs/>`_ and is tested on the following:

- Ubuntu 18, gcc 7.5.0, Python 3.6.9
- Ubuntu 18, gcc 7.5.0, Python 3.6.9 and 3.7
- Windows 10, VS 2019, Python 3.8.2 x64

It has dependencies on numpy and Cython that require these to be installed before installation of scikit-highs. Once numpy and Cython are installed, install using pip like this:
Expand All @@ -18,8 +18,8 @@ Example usage from script:

.. code:: python
from pyHiGHS import highs_wrapper
from pyHiGHS import CONST_INF, CONST_I_INF
from _highs.highs_wrapper import highs_wrapper
from _highs.highs_wrapper import CONST_INF, CONST_I_INF
# Call like this (usually using numpy arrays and A is a scipy.sparse.csc_matrix):
res = highs_wrapper(c, A.indptr, A.indices, A.data, lhs, rhs, lb, ub, options)
Expand Down
7 changes: 7 additions & 0 deletions scipy/optimize/_highs/pyHiGHS/src/HConst.pxd
Original file line number Diff line number Diff line change
Expand Up @@ -30,3 +30,10 @@ cdef extern from "HConst.h" nogil:
SOLVER_OPTION_SIMPLEX "SolverOption::SOLVER_OPTION_SIMPLEX" = -1
SOLVER_OPTION_CHOOSE "SolverOption::SOLVER_OPTION_CHOOSE"
SOLVER_OPTION_IPM "SolverOption::SOLVER_OPTION_IPM"

cdef enum PrimalDualStatus:
PrimalDualStatusSTATUS_NOT_SET "PrimalDualStatus::STATUS_NOT_SET" = -1
PrimalDualStatusSTATUS_NO_SOLUTION "PrimalDualStatus::STATUS_NO_SOLUTION"
PrimalDualStatusSTATUS_UNKNOWN "PrimalDualStatus::STATUS_UNKOWN"
PrimalDualStatusSTATUS_INFEASIBLE_POINT "PrimalDualStatus::STATUS_INFEASIBLE_POINT"
PrimalDualStatusSTATUS_FEASIBLE_POINT "PrimalDualStatus::STATUS_FEASIBLE_POINT"
4 changes: 4 additions & 0 deletions scipy/optimize/_highs/pyHiGHS/src/HighsLp.pxd
Original file line number Diff line number Diff line change
Expand Up @@ -40,6 +40,10 @@ cdef extern from "HighsLp.h" nogil:
HighsModelStatusREACHED_TIME_LIMIT "HighsModelStatus::REACHED_TIME_LIMIT"
HighsModelStatusREACHED_ITERATION_LIMIT "HighsModelStatus::REACHED_ITERATION_LIMIT"

ctypedef enum ObjSense:
ObjSenseMINIMIZE "ObjSense::MINIMIZE" = 1
ObjSenseMAXIMIZE "ObjSense::MAXIMIZE" = -1

cdef cppclass HighsSolution:
vector[double] col_value
vector[double] col_dual
Expand Down
19 changes: 8 additions & 11 deletions scipy/optimize/_highs/pyHiGHS/src/highs_wrapper.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -6,7 +6,10 @@ import logging
from libc.stdio cimport FILE, tmpfile
from libcpp.memory cimport unique_ptr

from HConst cimport ML_NONE
from HConst cimport (
ML_NONE,
PrimalDualStatusSTATUS_FEASIBLE_POINT,
)
from Highs cimport Highs
from HighsStatus cimport (
HighsStatus,
Expand Down Expand Up @@ -392,7 +395,7 @@ def highs_wrapper(
model_status = scaled_model_status

# We might need an info object if we can look up the solution and a place to put solution
cdef HighsInfo info
cdef HighsInfo info = highs.getHighsInfo() # it should always be safe to get the info object
cdef HighsSolution solution

# If the status is bad, don't look up the solution
Expand All @@ -408,23 +411,17 @@ def highs_wrapper(
HighsModelStatusPRIMAL_UNBOUNDED,
HighsModelStatusREACHED_ITERATION_LIMIT,
]:
info = highs.getHighsInfo()
return {
'status': <int> model_status,
'message': highs.highsModelStatusToString(model_status).decode(),
'simplex_nit': info.simplex_iteration_count,
'ipm_nit': info.ipm_iteration_count,
'fun': info.objective_function_value,
'crossover_nit': info.crossover_iteration_count,
'con': info.sum_primal_infeasibilities,
}
# If the model status is such that the solution can be read
elif model_status in [
HighsModelStatusOPTIMAL,
HighsModelStatusREACHED_DUAL_OBJECTIVE_VALUE_UPPER_BOUND,
HighsModelStatusREACHED_TIME_LIMIT,
]:
info = highs.getHighsInfo()
else:
# Should be safe to read the solution:
solution = highs.getSolution()
return {
'status': <int> model_status,
Expand All @@ -434,6 +431,7 @@ def highs_wrapper(
'x': [solution.col_value[ii] for ii in range(numcol)],

# Ax + s = b => Ax = b - s
# Note: this is for all constraints (A_ub and A_eq)
'slack': [rhs[ii] - solution.row_value[ii] for ii in range(numrow)],

# slacks in HiGHS appear as Ax - s, not Ax + s, so lambda is negated;
Expand All @@ -447,7 +445,6 @@ def highs_wrapper(
'simplex_nit': info.simplex_iteration_count,
'ipm_nit': info.ipm_iteration_count,
'crossover_nit': info.crossover_iteration_count,
'con': info.sum_primal_infeasibilities,
}

# Export some things
Expand Down
45 changes: 39 additions & 6 deletions scipy/optimize/_highs/pyHiGHS/src/mpswriter.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -7,12 +7,16 @@ from libcpp.vector cimport vector
from libcpp.string cimport string

from HighsStatus cimport HighsStatus
from HighsLp cimport (
ObjSense,
ObjSenseMINIMIZE,
)

from scipy.sparse import csc_matrix

cdef extern from "HMPSIO.h" nogil:
HighsStatus writeMPS(FILE* logfile, const char* filename, const int& numRow, const int& numCol,
const int& numInt, const int& objSense, const double& objOffset,
const int& numInt, const ObjSense& objSense, const double& objOffset,
const vector[int]& Astart, const vector[int]& Aindex,
const vector[double]& Avalue, const vector[double]& colCost,
const vector[double]& colLower, const vector[double]& colUpper,
Expand All @@ -24,20 +28,49 @@ def mpswriter(
string filename,
const double[::1] c,
A,
const double[::1] rhs,
const double[::1] lhs,
const double[::1] rhs,
const double[::1] lb,
const double[::1] ub,
const int[::1] integer_valued,
const bool use_free_format=True):
const bool use_free_format=False):
'''MPS writer: create MPS file from matrices.
Parameters
----------
filename : bytes (string)
Name of MPS file to write out to. Will be overwritten.
c : 1-D array (numcol,)
Objective coefficient values (assumes minimization).
A : 2-D array (numrow, numcol), scipy.sparse.csc_matrix
Sparse inequality constraint matrix.
lhs : 1-D array (numrow,)
Left hand side inequality values.
rhs : 1-D array (numrow,)
Right hand side inequality values.
lb : 1-D array (numcol,)
Lower bounds of solution variables.
ub : 1-D array (numcol,)
Upper bounds of solution variables.
integer_valued : 1-D array (numint,)
Indices of integer valued solution variables.
use_free_format : bool, optional
Use MPS free format. Default is False.
Notes
-----
Wrapper over HiGHS `writeMPS` function.
'''

if not isinstance(A, csc_matrix):
raise ValueError('A must be a scipy.sparse.csc_matrix!')

cdef int ii = 0
A = csc_matrix(A)

cdef int numRow = A.shape[0]
cdef int numCol = A.shape[1]
cdef int numInt = A.nnz
cdef int objSense = 1 # MIN for now
cdef int numInt = integer_valued.size
cdef ObjSense objSense = ObjSenseMINIMIZE # MIN for now
cdef double objOffset = 0 # This is a RHS on cost row

cdef vector[int] Astart = A.indptr
Expand Down
66 changes: 45 additions & 21 deletions scipy/optimize/_highs/setup.py
Original file line number Diff line number Diff line change
Expand Up @@ -25,7 +25,6 @@ def _get_sources(CMakeLists, start_token, end_token):
sources = [str(pathlib.Path('src/' + s)) for s in sources]
return sources


# Grab some more info about HiGHS from root CMakeLists
def _get_version(CMakeLists, start_token, end_token=')'):
CMakeLists = pathlib.Path(__file__).parent / CMakeLists
Expand Down Expand Up @@ -80,36 +79,35 @@ def configuration(parent_package='', top_path=None):
'basiclu',
sources=basiclu_sources,
include_dirs=[
str(pathlib.Path('src/')),
str(pathlib.Path('src/ipm/basiclu/include/')),
'src/',
'src/ipm/basiclu/include/',
],
language='c',
macros=DEFINE_MACROS,
)

# Compile the rest of the sources all together,
# linking the BASICLU static library
# highs_wrapper:
ipx_sources = _get_sources('src/CMakeLists.txt', 'set(ipx_sources\n', ')')
highs_sources = _get_sources('src/CMakeLists.txt', 'set(sources\n', ')')
WRAPPER_INCLUDE_DIRS = [
str(pathlib.Path('src/ipm/basiclu/include/')),
str(pathlib.Path('external/')),
str(pathlib.Path('src/')),
str(pathlib.Path('src/ipm/ipx/include/')),
str(pathlib.Path('src/lp_data/')),
str(pathlib.Path('src/io/')),
str(pathlib.Path('src/mip/')),
str(pathlib.Path('src/interfaces/')),
str(pathlib.Path('pyHiGHS/src/')),
]
ext = config.add_extension(
'highs_wrapper',
sources=[
str(pathlib.Path('pyHiGHS/src/highs_wrapper.cxx'))
] + ipx_sources + highs_sources,
sources=['pyHiGHS/src/highs_wrapper.cxx'] + highs_sources + ipx_sources,
include_dirs=[
str(pathlib.Path('pyHiGHS/src/')),
] + WRAPPER_INCLUDE_DIRS,

# highs_wrapper
'pyHiGHS/src/',
'src/',
'src/lp_data/',

# highs
'src/',
'src/io/',
'src/ipm/ipx/include/',

# IPX
'src/ipm/ipx/include/',
'src/ipm/basiclu/include/',
],
language='c++',
libraries=['basiclu'],
define_macros=DEFINE_MACROS,
Expand All @@ -118,6 +116,32 @@ def configuration(parent_package='', top_path=None):
# Add c++11/14 support:
ext._pre_build_hook = pre_build_hook

# wrapper around HiGHS writeMPS:
ext = config.add_extension(
'mpswriter',
sources=[
# we should be using using highs shared library;
# next best thing is compiling minimal set of sources
'pyHiGHS/src/mpswriter.cxx',
'src/util/HighsUtils.cpp',
'src/io/HighsIO.cpp',
'src/io/HMPSIO.cpp',
'src/lp_data/HighsModelUtils.cpp',
'src/util/stringutil.cpp',
],
include_dirs=[
'pyHiGHS/src/',
'src/',
'src/io/',
'src/lp_data/',
],
language='c++',
libraries=['basiclu'],
define_macros=DEFINE_MACROS,
undef_macros=UNDEF_MACROS,
)
ext._pre_build_hook = pre_build_hook

return config


Expand Down

0 comments on commit 2aa31fa

Please sign in to comment.