#include <sstream>
#include "Highs.h"
#include "lp_data/HighsLpUtils.h"
#include "lp_data/HighsModelUtils.h"
#include "mip/HighsMipSolver.h"
#include "model/HighsHessianUtils.h"
#include "simplex/HSimplex.h"
#include "util/HighsMatrixUtils.h"
#include "util/HighsSort.h"
void Highs::reportModelStats() const {
const HighsLp& lp = this->model_.lp_;
const HighsHessian& hessian = this->model_.hessian_;
const HighsLogOptions& log_options = this->options_.log_options;
if (!*log_options.output_flag) return;
HighsInt num_integer = 0;
HighsInt num_binary = 0;
HighsInt num_semi_continuous = 0;
HighsInt num_semi_integer = 0;
for (HighsInt iCol = 0; iCol < static_cast<HighsInt>(lp.integrality_.size());
iCol++) {
switch (lp.integrality_[iCol]) {
case HighsVarType::kInteger:
num_integer++;
if (lp.col_lower_[iCol] == 0 && lp.col_upper_[iCol] == 1) num_binary++;
break;
case HighsVarType::kSemiContinuous:
num_semi_continuous++;
break;
case HighsVarType::kSemiInteger:
num_semi_integer++;
break;
default:
break;
}
}
std::string problem_type;
const bool non_continuous =
num_integer + num_semi_continuous + num_semi_integer;
if (hessian.dim_) {
if (non_continuous) {
problem_type = "MIQP";
} else {
problem_type = "QP";
}
} else {
if (non_continuous) {
problem_type = "MIP";
} else {
problem_type = "LP";
}
}
const HighsInt a_num_nz = lp.a_matrix_.numNz();
const HighsInt q_num_nz = hessian.dim_ > 0 ? hessian.numNz() : 0;
if (*log_options.log_dev_level) {
highsLogDev(log_options, HighsLogType::kInfo, "%4s : %s\n",
problem_type.c_str(), lp.model_name_.c_str());
highsLogDev(log_options, HighsLogType::kInfo,
"Row%s : %" HIGHSINT_FORMAT "\n",
lp.num_row_ == 1 ? "" : "s", lp.num_row_);
highsLogDev(log_options, HighsLogType::kInfo,
"Col%s : %" HIGHSINT_FORMAT "\n",
lp.num_col_ == 1 ? "" : "s", lp.num_col_);
if (q_num_nz) {
highsLogDev(log_options, HighsLogType::kInfo,
"Matrix Nz : %" HIGHSINT_FORMAT "\n", a_num_nz);
highsLogDev(log_options, HighsLogType::kInfo,
"Hessian Nz: %" HIGHSINT_FORMAT "\n", q_num_nz);
} else {
highsLogDev(log_options, HighsLogType::kInfo,
"Nonzero%s : %" HIGHSINT_FORMAT "\n",
a_num_nz == 1 ? "" : "s", a_num_nz);
}
if (num_integer)
highsLogDev(log_options, HighsLogType::kInfo,
"Integer : %" HIGHSINT_FORMAT " (%" HIGHSINT_FORMAT
" binary)\n",
num_integer, num_binary);
if (num_semi_continuous)
highsLogDev(log_options, HighsLogType::kInfo,
"SemiConts : %" HIGHSINT_FORMAT "\n", num_semi_continuous);
if (num_semi_integer)
highsLogDev(log_options, HighsLogType::kInfo,
"SemiInt : %" HIGHSINT_FORMAT "\n", num_semi_integer);
} else {
std::stringstream stats_line;
stats_line << problem_type;
if (lp.model_name_.length()) stats_line << " " << lp.model_name_;
stats_line << " has " << lp.num_row_ << " row"
<< (lp.num_row_ == 1 ? "" : "s") << "; " << lp.num_col_ << " col"
<< (lp.num_col_ == 1 ? "" : "s");
if (q_num_nz) {
stats_line << "; " << a_num_nz << " matrix nonzero"
<< (a_num_nz == 1 ? "" : "s");
stats_line << "; " << q_num_nz << " Hessian nonzero"
<< (q_num_nz == 1 ? "" : "s");
} else {
stats_line << "; " << a_num_nz << " nonzero"
<< (a_num_nz == 1 ? "" : "s");
}
if (num_integer)
stats_line << "; " << num_integer << " integer variable"
<< (a_num_nz == 1 ? "" : "s") << " (" << num_binary
<< " binary)";
if (num_semi_continuous)
stats_line << "; " << num_semi_continuous << " semi-continuous variables";
if (num_semi_integer)
stats_line << "; " << num_semi_integer << " semi-integer variables";
highsLogUser(log_options, HighsLogType::kInfo, "%s\n",
stats_line.str().c_str());
}
}
HighsStatus Highs::formStandardFormLp() {
this->clearStandardFormLp();
HighsLp& lp = this->model_.lp_;
HighsSparseMatrix& matrix = lp.a_matrix_;
matrix.ensureRowwise();
HighsInt sense = HighsInt(lp.sense_);
this->standard_form_offset_ = sense * lp.offset_;
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++)
this->standard_form_cost_.push_back(sense * lp.col_cost_[iCol]);
this->standard_form_matrix_.format_ = MatrixFormat::kRowwise;
this->standard_form_matrix_.num_col_ = lp.num_col_;
HighsInt local_row_min_nnz = std::max(lp.num_col_, HighsInt(2));
HighsSparseMatrix local_row;
local_row.ensureRowwise();
local_row.num_row_ = 1;
local_row.num_col_ = lp.num_col_;
local_row.index_.resize(local_row_min_nnz);
local_row.value_.resize(local_row_min_nnz);
local_row.start_.resize(2);
HighsInt& num_nz = local_row.start_[1];
local_row.start_[0] = 0;
HighsInt num_fixed_row = 0;
HighsInt num_boxed_row = 0;
HighsInt num_lower_row = 0;
HighsInt num_upper_row = 0;
HighsInt num_free_row = 0;
HighsInt num_fixed_col = 0;
HighsInt num_boxed_col = 0;
HighsInt num_lower_col = 0;
HighsInt num_upper_col = 0;
HighsInt num_free_col = 0;
std::vector<HighsInt> slack_ix;
for (HighsInt iRow = 0; iRow < lp.num_row_; iRow++) {
double lower = lp.row_lower_[iRow];
double upper = lp.row_upper_[iRow];
if (lower <= -kHighsInf && upper >= kHighsInf) {
assert(0 == 1);
num_free_row++;
continue;
}
if (lower == upper) {
num_fixed_row++;
matrix.getRow(iRow, num_nz, local_row.index_.data(),
local_row.value_.data());
this->standard_form_matrix_.addRows(local_row);
this->standard_form_rhs_.push_back(upper);
continue;
} else if (lower <= -kHighsInf) {
num_upper_row++;
assert(upper < kHighsInf);
HighsInt standard_form_row = this->standard_form_rhs_.size();
slack_ix.push_back(standard_form_row + 1);
matrix.getRow(iRow, num_nz, local_row.index_.data(),
local_row.value_.data());
this->standard_form_matrix_.addRows(local_row);
this->standard_form_rhs_.push_back(upper);
} else if (upper >= kHighsInf) {
num_lower_row++;
assert(lower > -kHighsInf);
HighsInt standard_form_row = this->standard_form_rhs_.size();
slack_ix.push_back(-(standard_form_row + 1));
matrix.getRow(iRow, num_nz, local_row.index_.data(),
local_row.value_.data());
this->standard_form_matrix_.addRows(local_row);
this->standard_form_rhs_.push_back(lower);
} else {
assert(lower > -kHighsInf);
assert(upper < kHighsInf);
num_boxed_row++;
HighsInt standard_form_row = this->standard_form_rhs_.size();
slack_ix.push_back(-(standard_form_row + 1));
matrix.getRow(iRow, num_nz, local_row.index_.data(),
local_row.value_.data());
this->standard_form_matrix_.addRows(local_row);
this->standard_form_rhs_.push_back(lower);
standard_form_row = this->standard_form_rhs_.size();
slack_ix.push_back(standard_form_row + 1);
this->standard_form_matrix_.addRows(local_row);
this->standard_form_rhs_.push_back(upper);
}
}
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
double lower = lp.col_lower_[iCol];
double upper = lp.col_upper_[iCol];
if (lower > -kHighsInf && upper < kHighsInf) {
this->standard_form_cost_.push_back(0);
this->standard_form_matrix_.num_col_++;
local_row.num_col_++;
local_row.index_[0] = iCol;
local_row.index_[1] = this->standard_form_matrix_.num_col_ - 1;
local_row.value_[0] = 1;
local_row.value_[1] = 1;
local_row.start_[1] = 2;
this->standard_form_matrix_.addRows(local_row);
this->standard_form_rhs_.push_back(upper);
}
}
matrix.ensureColwise();
this->standard_form_matrix_.ensureColwise();
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
double cost = sense * lp.col_cost_[iCol];
double lower = lp.col_lower_[iCol];
double upper = lp.col_upper_[iCol];
if (lower > -kHighsInf) {
if (upper < kHighsInf) {
if (lower == upper) {
num_fixed_col++;
} else {
num_boxed_col++;
}
} else {
num_lower_col++;
}
if (lower != 0) {
this->standard_form_offset_ += cost * lower;
for (HighsInt iEl = this->standard_form_matrix_.start_[iCol];
iEl < this->standard_form_matrix_.start_[iCol + 1]; iEl++)
this->standard_form_rhs_[this->standard_form_matrix_.index_[iEl]] -=
this->standard_form_matrix_.value_[iEl] * lower;
}
} else if (upper < kHighsInf) {
num_upper_col++;
this->standard_form_offset_ += cost * upper;
this->standard_form_cost_[iCol] = -cost;
for (HighsInt iEl = this->standard_form_matrix_.start_[iCol];
iEl < this->standard_form_matrix_.start_[iCol + 1]; iEl++) {
this->standard_form_rhs_[this->standard_form_matrix_.index_[iEl]] -=
this->standard_form_matrix_.value_[iEl] * upper;
this->standard_form_matrix_.value_[iEl] =
-this->standard_form_matrix_.value_[iEl];
}
} else {
num_free_col++;
this->standard_form_cost_.push_back(-cost);
for (HighsInt iEl = this->standard_form_matrix_.start_[iCol];
iEl < this->standard_form_matrix_.start_[iCol + 1]; iEl++) {
this->standard_form_matrix_.index_.push_back(
this->standard_form_matrix_.index_[iEl]);
this->standard_form_matrix_.value_.push_back(
-this->standard_form_matrix_.value_[iEl]);
}
this->standard_form_matrix_.start_.push_back(
HighsInt(this->standard_form_matrix_.index_.size()));
}
}
for (HighsInt iRow : slack_ix) {
this->standard_form_cost_.push_back(0);
if (iRow > 0) {
this->standard_form_matrix_.index_.push_back(iRow - 1);
this->standard_form_matrix_.value_.push_back(1);
} else {
this->standard_form_matrix_.index_.push_back(-iRow - 1);
this->standard_form_matrix_.value_.push_back(-1);
}
this->standard_form_matrix_.start_.push_back(
HighsInt(this->standard_form_matrix_.index_.size()));
}
this->standard_form_matrix_.num_col_ = int(standard_form_cost_.size());
this->standard_form_matrix_.num_row_ = int(standard_form_rhs_.size());
this->standard_form_valid_ = true;
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Standard form LP obtained for LP with (free / lower / upper / "
"boxed / fixed) variables"
" (%d / %d / %d / %d / %d) and constraints"
" (%d / %d / %d / %d / %d) \n",
int(num_free_col), int(num_lower_col), int(num_upper_col),
int(num_boxed_col), int(num_fixed_col), int(num_free_row),
int(num_lower_row), int(num_upper_row), int(num_boxed_row),
int(num_fixed_row));
return HighsStatus::kOk;
}
HighsStatus Highs::basisForSolution() {
HighsLp& lp = model_.lp_;
assert(!lp.isMip() || options_.solve_relaxation);
assert(solution_.value_valid);
invalidateBasis();
HighsInt num_basic = 0;
HighsBasis basis;
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
if (std::fabs(lp.col_lower_[iCol] - solution_.col_value[iCol]) <=
options_.primal_feasibility_tolerance) {
basis.col_status.push_back(HighsBasisStatus::kLower);
} else if (std::fabs(lp.col_upper_[iCol] - solution_.col_value[iCol]) <=
options_.primal_feasibility_tolerance) {
basis.col_status.push_back(HighsBasisStatus::kUpper);
} else {
num_basic++;
basis.col_status.push_back(HighsBasisStatus::kBasic);
}
}
const HighsInt num_basic_col = num_basic;
for (HighsInt iRow = 0; iRow < lp.num_row_; iRow++) {
if (std::fabs(lp.row_lower_[iRow] - solution_.row_value[iRow]) <=
options_.primal_feasibility_tolerance) {
basis.row_status.push_back(HighsBasisStatus::kLower);
} else if (std::fabs(lp.row_upper_[iRow] - solution_.row_value[iRow]) <=
options_.primal_feasibility_tolerance) {
basis.row_status.push_back(HighsBasisStatus::kUpper);
} else {
num_basic++;
basis.row_status.push_back(HighsBasisStatus::kBasic);
}
}
const HighsInt num_basic_row = num_basic - num_basic_col;
assert((int)basis.col_status.size() == lp.num_col_);
assert((int)basis.row_status.size() == lp.num_row_);
highsLogDev(options_.log_options, HighsLogType::kInfo,
"LP has %d rows and solution yields %d possible basic variables "
"(%d / %d; %d / %d)\n",
(int)lp.num_row_, (int)num_basic, (int)num_basic_col,
(int)lp.num_col_, (int)num_basic_row, (int)lp.num_row_);
return this->setBasis(basis);
}
HighsStatus Highs::addColsInterface(
HighsInt ext_num_new_col, const double* ext_col_cost,
const double* ext_col_lower, const double* ext_col_upper,
HighsInt ext_num_new_nz, const HighsInt* ext_a_start,
const HighsInt* ext_a_index, const double* ext_a_value) {
HighsStatus return_status = HighsStatus::kOk;
HighsOptions& options = options_;
if (ext_num_new_col < 0) return HighsStatus::kError;
if (ext_num_new_nz < 0) return HighsStatus::kError;
if (ext_num_new_col == 0) return HighsStatus::kOk;
if (ext_num_new_col > 0)
if (isColDataNull(options.log_options, ext_col_cost, ext_col_lower,
ext_col_upper))
return HighsStatus::kError;
if (ext_num_new_nz > 0)
if (isMatrixDataNull(options.log_options, ext_a_start, ext_a_index,
ext_a_value))
return HighsStatus::kError;
HighsLp& lp = model_.lp_;
HighsBasis& basis = basis_;
HighsScale& scale = lp.scale_;
bool& useful_basis = basis.useful;
bool& lp_has_scaling = lp.scale_.has_scaling;
if (lp.num_row_ <= 0 && ext_num_new_nz > 0) return HighsStatus::kError;
HighsInt newNumCol = lp.num_col_ + ext_num_new_col;
HighsIndexCollection index_collection;
index_collection.dimension_ = ext_num_new_col;
index_collection.is_interval_ = true;
index_collection.from_ = 0;
index_collection.to_ = ext_num_new_col - 1;
std::vector<double> local_colCost{ext_col_cost,
ext_col_cost + ext_num_new_col};
std::vector<double> local_colLower{ext_col_lower,
ext_col_lower + ext_num_new_col};
std::vector<double> local_colUpper{ext_col_upper,
ext_col_upper + ext_num_new_col};
bool local_has_infinite_cost = false;
return_status = interpretCallStatus(
options_.log_options,
assessCosts(options, lp.num_col_, index_collection, local_colCost,
local_has_infinite_cost, options.infinite_cost),
return_status, "assessCosts");
if (return_status == HighsStatus::kError) return return_status;
return_status = interpretCallStatus(
options_.log_options,
assessBounds(options, "Col", lp.num_col_, index_collection,
local_colLower, local_colUpper, options.infinite_bound),
return_status, "assessBounds");
if (return_status == HighsStatus::kError) return return_status;
appendColsToLpVectors(lp, ext_num_new_col, local_colCost, local_colLower,
local_colUpper);
HighsSparseMatrix local_a_matrix;
local_a_matrix.num_col_ = ext_num_new_col;
local_a_matrix.num_row_ = lp.num_row_;
local_a_matrix.format_ = MatrixFormat::kColwise;
if (ext_num_new_nz) {
local_a_matrix.start_ = {ext_a_start, ext_a_start + ext_num_new_col};
local_a_matrix.start_.resize(ext_num_new_col + 1);
local_a_matrix.start_[ext_num_new_col] = ext_num_new_nz;
local_a_matrix.index_ = {ext_a_index, ext_a_index + ext_num_new_nz};
local_a_matrix.value_ = {ext_a_value, ext_a_value + ext_num_new_nz};
return_status =
interpretCallStatus(options_.log_options,
local_a_matrix.assess(options.log_options, "LP",
options.small_matrix_value,
options.large_matrix_value),
return_status, "assessMatrix");
if (return_status == HighsStatus::kError) return return_status;
} else {
local_a_matrix.start_.assign(ext_num_new_col + 1, 0);
}
lp.a_matrix_.addCols(local_a_matrix);
if (lp_has_scaling) {
scale.col.resize(newNumCol);
for (HighsInt iCol = 0; iCol < ext_num_new_col; iCol++)
scale.col[lp.num_col_ + iCol] = 1.0;
scale.num_col = newNumCol;
local_a_matrix.applyRowScale(scale);
local_a_matrix.considerColScaling(options.allowed_matrix_scale_factor,
&scale.col[lp.num_col_]);
}
if (useful_basis) appendNonbasicColsToBasisInterface(ext_num_new_col);
lp.addColNames("", ext_num_new_col);
lp.num_col_ += ext_num_new_col;
assert(lpDimensionsOk("addCols", lp, options.log_options));
lp.has_infinite_cost_ = lp.has_infinite_cost_ || local_has_infinite_cost;
assert(lp.has_infinite_cost_ == lp.hasInfiniteCost(options.infinite_cost));
invalidateModelStatusSolutionAndInfo();
ekk_instance_.addCols(lp, local_a_matrix);
if (this->model_.hessian_.dim_)
completeHessian(lp.num_col_, this->model_.hessian_);
return return_status;
}
HighsStatus Highs::addRowsInterface(HighsInt ext_num_new_row,
const double* ext_row_lower,
const double* ext_row_upper,
HighsInt ext_num_new_nz,
const HighsInt* ext_ar_start,
const HighsInt* ext_ar_index,
const double* ext_ar_value) {
if (kExtendInvertWhenAddingRows) {
if (ekk_instance_.status_.has_nla)
ekk_instance_.debugNlaCheckInvert("Start of Highs::addRowsInterface",
kHighsDebugLevelExpensive + 1);
}
HighsStatus return_status = HighsStatus::kOk;
HighsOptions& options = options_;
if (ext_num_new_row < 0) return HighsStatus::kError;
if (ext_num_new_nz < 0) return HighsStatus::kError;
if (ext_num_new_row == 0) return HighsStatus::kOk;
if (ext_num_new_row > 0)
if (isRowDataNull(options.log_options, ext_row_lower, ext_row_upper))
return HighsStatus::kError;
if (ext_num_new_nz > 0)
if (isMatrixDataNull(options.log_options, ext_ar_start, ext_ar_index,
ext_ar_value))
return HighsStatus::kError;
HighsLp& lp = model_.lp_;
HighsBasis& basis = basis_;
HighsScale& scale = lp.scale_;
bool& useful_basis = basis.useful;
bool& lp_has_scaling = lp.scale_.has_scaling;
if (lp.num_col_ <= 0 && ext_num_new_nz > 0) return HighsStatus::kError;
HighsInt newNumRow = lp.num_row_ + ext_num_new_row;
HighsIndexCollection index_collection;
index_collection.dimension_ = ext_num_new_row;
index_collection.is_interval_ = true;
index_collection.from_ = 0;
index_collection.to_ = ext_num_new_row - 1;
std::vector<double> local_rowLower{ext_row_lower,
ext_row_lower + ext_num_new_row};
std::vector<double> local_rowUpper{ext_row_upper,
ext_row_upper + ext_num_new_row};
return_status = interpretCallStatus(
options_.log_options,
assessBounds(options, "Row", lp.num_row_, index_collection,
local_rowLower, local_rowUpper, options.infinite_bound),
return_status, "assessBounds");
if (return_status == HighsStatus::kError) return return_status;
appendRowsToLpVectors(lp, ext_num_new_row, local_rowLower, local_rowUpper);
HighsSparseMatrix local_ar_matrix;
local_ar_matrix.num_col_ = lp.num_col_;
local_ar_matrix.num_row_ = ext_num_new_row;
local_ar_matrix.format_ = MatrixFormat::kRowwise;
if (ext_num_new_nz) {
local_ar_matrix.start_ = {ext_ar_start, ext_ar_start + ext_num_new_row};
local_ar_matrix.start_.resize(ext_num_new_row + 1);
local_ar_matrix.start_[ext_num_new_row] = ext_num_new_nz;
local_ar_matrix.index_ = {ext_ar_index, ext_ar_index + ext_num_new_nz};
local_ar_matrix.value_ = {ext_ar_value, ext_ar_value + ext_num_new_nz};
return_status =
interpretCallStatus(options_.log_options,
local_ar_matrix.assess(options.log_options, "LP",
options.small_matrix_value,
options.large_matrix_value),
return_status, "assessMatrix");
if (return_status == HighsStatus::kError) return return_status;
} else {
local_ar_matrix.start_.assign(ext_num_new_row + 1, 0);
}
lp.a_matrix_.addRows(local_ar_matrix);
if (lp_has_scaling) {
scale.row.resize(newNumRow);
for (HighsInt iRow = 0; iRow < ext_num_new_row; iRow++)
scale.row[lp.num_row_ + iRow] = 1.0;
scale.num_row = newNumRow;
local_ar_matrix.applyColScale(scale);
local_ar_matrix.considerRowScaling(options.allowed_matrix_scale_factor,
&scale.row[lp.num_row_]);
}
if (useful_basis) appendBasicRowsToBasisInterface(ext_num_new_row);
lp.addRowNames("", ext_num_new_row);
lp.num_row_ += ext_num_new_row;
assert(lpDimensionsOk("addRows", lp, options.log_options));
invalidateModelStatusSolutionAndInfo();
ekk_instance_.addRows(lp, local_ar_matrix);
return return_status;
}
static void deleteBasisEntries(std::vector<HighsBasisStatus>& status,
bool& deleted_basic, bool& deleted_nonbasic,
const HighsIndexCollection& index_collection,
const HighsInt entry_dim) {
assert(ok(index_collection));
assert(static_cast<size_t>(entry_dim) == status.size());
HighsInt from_k;
HighsInt to_k;
limits(index_collection, from_k, to_k);
if (from_k > to_k) return;
HighsInt delete_from_entry;
HighsInt delete_to_entry;
HighsInt keep_from_entry;
HighsInt keep_to_entry = -1;
HighsInt current_set_entry = 0;
HighsInt new_num_entry = 0;
deleted_basic = false;
deleted_nonbasic = false;
for (HighsInt k = from_k; k <= to_k; k++) {
updateOutInIndex(index_collection, delete_from_entry, delete_to_entry,
keep_from_entry, keep_to_entry, current_set_entry);
if (k == from_k) new_num_entry = delete_from_entry;
for (HighsInt entry = delete_from_entry; entry <= delete_to_entry;
entry++) {
if (status[entry] == HighsBasisStatus::kBasic) {
deleted_basic = true;
} else {
deleted_nonbasic = true;
}
}
if (delete_to_entry >= entry_dim - 1) break;
for (HighsInt entry = keep_from_entry; entry <= keep_to_entry; entry++) {
status[new_num_entry] = status[entry];
new_num_entry++;
}
if (keep_to_entry >= entry_dim - 1) break;
}
status.resize(new_num_entry);
}
static void deleteBasisCols(HighsBasis& basis,
const HighsIndexCollection& index_collection,
const HighsInt original_num_col) {
bool deleted_basic;
bool deleted_nonbasic;
deleteBasisEntries(basis.col_status, deleted_basic, deleted_nonbasic,
index_collection, original_num_col);
if (deleted_basic) basis.valid = false;
}
static void deleteBasisRows(HighsBasis& basis,
const HighsIndexCollection& index_collection,
const HighsInt original_num_row) {
bool deleted_basic;
bool deleted_nonbasic;
deleteBasisEntries(basis.row_status, deleted_basic, deleted_nonbasic,
index_collection, original_num_row);
if (deleted_nonbasic) basis.valid = false;
}
void Highs::deleteColsInterface(HighsIndexCollection& index_collection) {
HighsLp& lp = model_.lp_;
HighsBasis& basis = basis_;
lp.ensureColwise();
HighsInt original_num_col = lp.num_col_;
lp.deleteCols(index_collection);
model_.hessian_.deleteCols(index_collection);
if (lp.num_col_ == original_num_col) return;
assert(lp.num_col_ < original_num_col);
model_status_ = HighsModelStatus::kNotset;
if (basis_.useful) {
assert(basis_.col_status.size() == static_cast<size_t>(original_num_col));
deleteBasisCols(basis_, index_collection, original_num_col);
} else {
assert(!basis.valid);
}
if (lp.scale_.has_scaling) {
deleteScale(lp.scale_.col, index_collection);
lp.scale_.col.resize(lp.num_col_);
lp.scale_.num_col = lp.num_col_;
}
invalidateModelStatusSolutionAndInfo();
ekk_instance_.deleteCols(index_collection);
if (index_collection.is_mask_) {
HighsInt new_col = 0;
for (HighsInt col = 0; col < original_num_col; col++) {
if (!index_collection.mask_[col]) {
index_collection.mask_[col] = new_col;
new_col++;
} else {
index_collection.mask_[col] = -1;
}
}
assert(new_col == lp.num_col_);
}
assert(lpDimensionsOk("deleteCols", lp, options_.log_options));
lp.col_hash_.name2index.clear();
}
void Highs::deleteRowsInterface(HighsIndexCollection& index_collection) {
HighsLp& lp = model_.lp_;
HighsBasis& basis = basis_;
lp.ensureColwise();
HighsInt original_num_row = lp.num_row_;
lp.deleteRows(index_collection);
if (lp.num_row_ == original_num_row) return;
assert(lp.num_row_ < original_num_row);
model_status_ = HighsModelStatus::kNotset;
if (basis_.useful) {
assert(basis_.row_status.size() == static_cast<size_t>(original_num_row));
deleteBasisRows(basis_, index_collection, original_num_row);
} else {
assert(!basis.valid);
}
if (lp.scale_.has_scaling) {
deleteScale(lp.scale_.row, index_collection);
lp.scale_.row.resize(lp.num_row_);
lp.scale_.num_row = lp.num_row_;
}
invalidateModelStatusSolutionAndInfo();
ekk_instance_.deleteRows(index_collection);
if (index_collection.is_mask_) {
HighsInt new_row = 0;
for (HighsInt row = 0; row < original_num_row; row++) {
if (!index_collection.mask_[row]) {
index_collection.mask_[row] = new_row;
new_row++;
} else {
index_collection.mask_[row] = -1;
}
}
assert(new_row == lp.num_row_);
}
assert(lpDimensionsOk("deleteRows", lp, options_.log_options));
lp.row_hash_.name2index.clear();
}
void Highs::getColsInterface(const HighsIndexCollection& index_collection,
HighsInt& num_col, double* cost, double* lower,
double* upper, HighsInt& num_nz, HighsInt* start,
HighsInt* index, double* value) const {
const HighsLp& lp = model_.lp_;
if (lp.a_matrix_.isColwise()) {
getSubVectors(index_collection, lp.num_col_, lp.col_cost_.data(),
lp.col_lower_.data(), lp.col_upper_.data(), lp.a_matrix_,
num_col, cost, lower, upper, num_nz, start, index, value);
} else {
getSubVectorsTranspose(index_collection, lp.num_col_, lp.col_cost_.data(),
lp.col_lower_.data(), lp.col_upper_.data(),
lp.a_matrix_, num_col, cost, lower, upper, num_nz,
start, index, value);
}
}
void Highs::getRowsInterface(const HighsIndexCollection& index_collection,
HighsInt& num_row, double* lower, double* upper,
HighsInt& num_nz, HighsInt* start, HighsInt* index,
double* value) const {
const HighsLp& lp = model_.lp_;
if (lp.a_matrix_.isColwise()) {
getSubVectorsTranspose(index_collection, lp.num_row_, nullptr,
lp.row_lower_.data(), lp.row_upper_.data(),
lp.a_matrix_, num_row, nullptr, lower, upper, num_nz,
start, index, value);
} else {
getSubVectors(index_collection, lp.num_row_, nullptr, lp.row_lower_.data(),
lp.row_upper_.data(), lp.a_matrix_, num_row, nullptr, lower,
upper, num_nz, start, index, value);
}
}
void Highs::getCoefficientInterface(const HighsInt ext_row,
const HighsInt ext_col,
double& value) const {
const HighsLp& lp = model_.lp_;
assert(0 <= ext_row && ext_row < lp.num_row_);
assert(0 <= ext_col && ext_col < lp.num_col_);
value = 0;
if (lp.a_matrix_.isColwise()) {
for (HighsInt el = lp.a_matrix_.start_[ext_col];
el < lp.a_matrix_.start_[ext_col + 1]; el++) {
if (lp.a_matrix_.index_[el] == ext_row) {
value = lp.a_matrix_.value_[el];
break;
}
}
} else {
for (HighsInt el = lp.a_matrix_.start_[ext_row];
el < lp.a_matrix_.start_[ext_row + 1]; el++) {
if (lp.a_matrix_.index_[el] == ext_col) {
value = lp.a_matrix_.value_[el];
break;
}
}
}
}
HighsStatus Highs::changeIntegralityInterface(
HighsIndexCollection& index_collection, const HighsVarType* integrality) {
HighsInt num_integrality = dataSize(index_collection);
if (num_integrality <= 0) return HighsStatus::kOk;
if (highsVarTypeUserDataNotNull(options_.log_options, integrality,
"column integrality"))
return HighsStatus::kError;
std::vector<HighsVarType> local_integrality{integrality,
integrality + num_integrality};
if (index_collection.is_set_)
assert(increasingSetOk(index_collection.set_, 0,
index_collection.dimension_, true));
HighsStatus return_status = changeLpIntegrality(model_.lp_, index_collection,
local_integrality, options_);
assert(return_status != HighsStatus::kError);
invalidateModelStatus();
return return_status;
}
HighsStatus Highs::changeCostsInterface(HighsIndexCollection& index_collection,
const double* cost) {
HighsInt num_cost = dataSize(index_collection);
if (num_cost <= 0) return HighsStatus::kOk;
if (doubleUserDataNotNull(options_.log_options, cost, "column costs"))
return HighsStatus::kError;
std::vector<double> local_colCost{cost, cost + num_cost};
HighsStatus return_status = HighsStatus::kOk;
bool local_has_infinite_cost = false;
return_status = interpretCallStatus(
options_.log_options,
assessCosts(options_, 0, index_collection, local_colCost,
local_has_infinite_cost, options_.infinite_cost),
return_status, "assessCosts");
if (return_status == HighsStatus::kError) return return_status;
HighsLp& lp = model_.lp_;
changeLpCosts(lp, index_collection, local_colCost, options_.infinite_cost);
lp.has_infinite_cost_ = lp.has_infinite_cost_ || local_has_infinite_cost;
assert(lp.has_infinite_cost_ == lp.hasInfiniteCost(options_.infinite_cost));
invalidateModelStatusSolutionAndInfo();
ekk_instance_.updateStatus(LpAction::kNewCosts);
return HighsStatus::kOk;
}
bool Highs::feasibleWrtBounds(const bool columns) const {
if (this->info_.primal_solution_status != kSolutionStatusFeasible)
return false;
const HighsLp& lp = model_.lp_;
const double primal_feasibility_tolerance =
this->options_.primal_feasibility_tolerance;
std::vector<double> value =
columns ? this->solution_.col_value : this->solution_.row_value;
std::vector<double> lower = columns ? lp.col_lower_ : lp.row_lower_;
std::vector<double> upper = columns ? lp.col_upper_ : lp.row_upper_;
HighsInt dim = columns ? lp.num_col_ : lp.num_row_;
for (HighsInt iX = 0; iX < dim; iX++) {
if (value[iX] < lower[iX] - primal_feasibility_tolerance) return false;
if (value[iX] > upper[iX] + primal_feasibility_tolerance) return false;
}
return true;
}
HighsStatus Highs::changeColBoundsInterface(
HighsIndexCollection& index_collection, const double* col_lower,
const double* col_upper) {
HighsInt num_col_bounds = dataSize(index_collection);
if (num_col_bounds <= 0) return HighsStatus::kOk;
bool null_data = false;
null_data = doubleUserDataNotNull(options_.log_options, col_lower,
"column lower bounds") ||
null_data;
null_data = doubleUserDataNotNull(options_.log_options, col_upper,
"column upper bounds") ||
null_data;
if (null_data) return HighsStatus::kError;
std::vector<double> local_colLower{col_lower, col_lower + num_col_bounds};
std::vector<double> local_colUpper{col_upper, col_upper + num_col_bounds};
if (index_collection.is_set_)
sortSetData(index_collection.set_num_entries_, index_collection.set_,
col_lower, col_upper, NULL, local_colLower.data(),
local_colUpper.data(), NULL);
HighsStatus return_status = HighsStatus::kOk;
return_status = interpretCallStatus(
options_.log_options,
assessBounds(options_, "col", 0, index_collection, local_colLower,
local_colUpper, options_.infinite_bound),
return_status, "assessBounds");
if (return_status == HighsStatus::kError) return return_status;
HighsLp& lp = model_.lp_;
changeLpColBounds(lp, index_collection, local_colLower, local_colUpper);
setNonbasicStatusInterface(index_collection, true);
if (!this->basis_.useful && feasibleWrtBounds()) {
invalidateModelStatusAndInfo();
} else {
invalidateModelStatusSolutionAndInfo();
}
ekk_instance_.updateStatus(LpAction::kNewBounds);
return HighsStatus::kOk;
}
HighsStatus Highs::changeRowBoundsInterface(
HighsIndexCollection& index_collection, const double* lower,
const double* upper) {
HighsInt num_row_bounds = dataSize(index_collection);
if (num_row_bounds <= 0) return HighsStatus::kOk;
bool null_data = false;
null_data =
doubleUserDataNotNull(options_.log_options, lower, "row lower bounds") ||
null_data;
null_data =
doubleUserDataNotNull(options_.log_options, upper, "row upper bounds") ||
null_data;
if (null_data) return HighsStatus::kError;
std::vector<double> local_rowLower{lower, lower + num_row_bounds};
std::vector<double> local_rowUpper{upper, upper + num_row_bounds};
if (index_collection.is_set_)
sortSetData(index_collection.set_num_entries_, index_collection.set_, lower,
upper, NULL, local_rowLower.data(), local_rowUpper.data(),
NULL);
HighsStatus return_status = HighsStatus::kOk;
return_status = interpretCallStatus(
options_.log_options,
assessBounds(options_, "row", 0, index_collection, local_rowLower,
local_rowUpper, options_.infinite_bound),
return_status, "assessBounds");
if (return_status == HighsStatus::kError) return return_status;
HighsLp& lp = model_.lp_;
changeLpRowBounds(lp, index_collection, local_rowLower, local_rowUpper);
setNonbasicStatusInterface(index_collection, false);
if (!this->basis_.useful && feasibleWrtBounds(false)) {
invalidateModelStatusAndInfo();
} else {
invalidateModelStatusSolutionAndInfo();
}
ekk_instance_.updateStatus(LpAction::kNewBounds);
return HighsStatus::kOk;
}
void Highs::changeCoefficientInterface(const HighsInt ext_row,
const HighsInt ext_col,
const double ext_new_value) {
HighsLp& lp = model_.lp_;
lp.ensureColwise();
assert(0 <= ext_row && ext_row < lp.num_row_);
assert(0 <= ext_col && ext_col < lp.num_col_);
const bool zero_new_value =
std::fabs(ext_new_value) <= options_.small_matrix_value;
changeLpMatrixCoefficient(lp, ext_row, ext_col, ext_new_value,
zero_new_value);
const bool basic_column =
this->basis_.col_status[ext_col] == HighsBasisStatus::kBasic;
invalidateModelStatusSolutionAndInfo();
if (basic_column) {
this->basis_.was_alien = true;
this->basis_.alien = true;
}
ekk_instance_.updateStatus(LpAction::kNewRows);
}
HighsStatus Highs::scaleColInterface(const HighsInt col,
const double scale_value) {
HighsStatus return_status = HighsStatus::kOk;
HighsLp& lp = model_.lp_;
HighsBasis& basis = basis_;
HighsSimplexStatus& simplex_status = ekk_instance_.status_;
lp.ensureColwise();
if (col < 0) return HighsStatus::kError;
if (col >= lp.num_col_) return HighsStatus::kError;
if (!scale_value) return HighsStatus::kError;
return_status = interpretCallStatus(options_.log_options,
applyScalingToLpCol(lp, col, scale_value),
return_status, "applyScalingToLpCol");
if (return_status == HighsStatus::kError) return return_status;
if (scale_value < 0 && basis.valid) {
if (basis.col_status[col] == HighsBasisStatus::kLower) {
basis.col_status[col] = HighsBasisStatus::kUpper;
} else if (basis.col_status[col] == HighsBasisStatus::kUpper) {
basis.col_status[col] = HighsBasisStatus::kLower;
}
}
if (simplex_status.initialised_for_solve) {
SimplexBasis& simplex_basis = ekk_instance_.basis_;
if (scale_value < 0 && simplex_status.has_basis) {
if (simplex_basis.nonbasicMove_[col] == kNonbasicMoveUp) {
simplex_basis.nonbasicMove_[col] = kNonbasicMoveDn;
} else if (simplex_basis.nonbasicMove_[col] == kNonbasicMoveDn) {
simplex_basis.nonbasicMove_[col] = kNonbasicMoveUp;
}
}
}
invalidateModelStatusSolutionAndInfo();
ekk_instance_.updateStatus(LpAction::kScaledCol);
return HighsStatus::kOk;
}
HighsStatus Highs::scaleRowInterface(const HighsInt row,
const double scale_value) {
HighsStatus return_status = HighsStatus::kOk;
HighsLp& lp = model_.lp_;
HighsBasis& basis = basis_;
HighsSimplexStatus& simplex_status = ekk_instance_.status_;
lp.ensureColwise();
if (row < 0) return HighsStatus::kError;
if (row >= lp.num_row_) return HighsStatus::kError;
if (!scale_value) return HighsStatus::kError;
return_status = interpretCallStatus(options_.log_options,
applyScalingToLpRow(lp, row, scale_value),
return_status, "applyScalingToLpRow");
if (return_status == HighsStatus::kError) return return_status;
if (scale_value < 0 && basis.valid) {
if (basis.row_status[row] == HighsBasisStatus::kLower) {
basis.row_status[row] = HighsBasisStatus::kUpper;
} else if (basis.row_status[row] == HighsBasisStatus::kUpper) {
basis.row_status[row] = HighsBasisStatus::kLower;
}
}
if (simplex_status.initialised_for_solve) {
SimplexBasis& simplex_basis = ekk_instance_.basis_;
if (scale_value < 0 && simplex_status.has_basis) {
const HighsInt var = lp.num_col_ + row;
if (simplex_basis.nonbasicMove_[var] == kNonbasicMoveUp) {
simplex_basis.nonbasicMove_[var] = kNonbasicMoveDn;
} else if (simplex_basis.nonbasicMove_[var] == kNonbasicMoveDn) {
simplex_basis.nonbasicMove_[var] = kNonbasicMoveUp;
}
}
}
invalidateModelStatusSolutionAndInfo();
ekk_instance_.updateStatus(LpAction::kScaledRow);
return HighsStatus::kOk;
}
void Highs::setNonbasicStatusInterface(
const HighsIndexCollection& index_collection, const bool columns) {
HighsBasis& highs_basis = basis_;
if (!highs_basis.valid) return;
const bool has_simplex_basis = ekk_instance_.status_.has_basis;
SimplexBasis& simplex_basis = ekk_instance_.basis_;
HighsLp& lp = model_.lp_;
assert(ok(index_collection));
HighsInt from_k;
HighsInt to_k;
limits(index_collection, from_k, to_k);
HighsInt ix_dim;
if (columns) {
ix_dim = lp.num_col_;
} else {
ix_dim = lp.num_row_;
}
assert(0 <= from_k && to_k < ix_dim);
assert(from_k <= to_k);
HighsInt set_from_ix;
HighsInt set_to_ix;
HighsInt ignore_from_ix;
HighsInt ignore_to_ix = -1;
HighsInt current_set_entry = 0;
for (HighsInt k = from_k; k <= to_k; k++) {
updateOutInIndex(index_collection, set_from_ix, set_to_ix, ignore_from_ix,
ignore_to_ix, current_set_entry);
assert(set_to_ix < ix_dim);
assert(ignore_to_ix < ix_dim);
if (columns) {
for (HighsInt iCol = set_from_ix; iCol <= set_to_ix; iCol++) {
if (highs_basis.col_status[iCol] == HighsBasisStatus::kBasic) continue;
double lower = lp.col_lower_[iCol];
double upper = lp.col_upper_[iCol];
HighsBasisStatus status = highs_basis.col_status[iCol];
int8_t move = kIllegalMoveValue;
if (lower == upper) {
if (status == HighsBasisStatus::kNonbasic)
status = HighsBasisStatus::kLower;
move = kNonbasicMoveZe;
} else if (!highs_isInfinity(-lower)) {
if (!highs_isInfinity(upper)) {
if (status == HighsBasisStatus::kNonbasic) {
if (fabs(lower) < fabs(upper)) {
status = HighsBasisStatus::kLower;
move = kNonbasicMoveUp;
} else {
status = HighsBasisStatus::kUpper;
move = kNonbasicMoveDn;
}
} else if (status == HighsBasisStatus::kLower) {
move = kNonbasicMoveUp;
} else {
move = kNonbasicMoveDn;
}
} else {
status = HighsBasisStatus::kLower;
move = kNonbasicMoveUp;
}
} else if (!highs_isInfinity(upper)) {
status = HighsBasisStatus::kUpper;
move = kNonbasicMoveDn;
} else {
status = HighsBasisStatus::kZero;
move = kNonbasicMoveZe;
}
highs_basis.col_status[iCol] = status;
if (has_simplex_basis) {
assert(move != kIllegalMoveValue);
simplex_basis.nonbasicFlag_[iCol] = kNonbasicFlagTrue;
simplex_basis.nonbasicMove_[iCol] = move;
}
}
} else {
for (HighsInt iRow = set_from_ix; iRow <= set_to_ix; iRow++) {
if (highs_basis.row_status[iRow] == HighsBasisStatus::kBasic) continue;
double lower = lp.row_lower_[iRow];
double upper = lp.row_upper_[iRow];
HighsBasisStatus status = highs_basis.row_status[iRow];
int8_t move = kIllegalMoveValue;
if (lower == upper) {
if (status == HighsBasisStatus::kNonbasic)
status = HighsBasisStatus::kLower;
move = kNonbasicMoveZe;
} else if (!highs_isInfinity(-lower)) {
if (!highs_isInfinity(upper)) {
if (status == HighsBasisStatus::kNonbasic) {
if (fabs(lower) < fabs(upper)) {
status = HighsBasisStatus::kLower;
move = kNonbasicMoveDn;
} else {
status = HighsBasisStatus::kUpper;
move = kNonbasicMoveUp;
}
} else if (status == HighsBasisStatus::kLower) {
move = kNonbasicMoveDn;
} else {
move = kNonbasicMoveUp;
}
} else {
status = HighsBasisStatus::kLower;
move = kNonbasicMoveDn;
}
} else if (!highs_isInfinity(upper)) {
status = HighsBasisStatus::kUpper;
move = kNonbasicMoveUp;
} else {
status = HighsBasisStatus::kZero;
move = kNonbasicMoveZe;
}
highs_basis.row_status[iRow] = status;
if (has_simplex_basis) {
assert(move != kIllegalMoveValue);
simplex_basis.nonbasicFlag_[lp.num_col_ + iRow] = kNonbasicFlagTrue;
simplex_basis.nonbasicMove_[lp.num_col_ + iRow] = move;
}
}
}
if (ignore_to_ix >= ix_dim - 1) break;
}
}
void Highs::appendNonbasicColsToBasisInterface(const HighsInt ext_num_new_col) {
if (ext_num_new_col == 0) return;
HighsBasis& highs_basis = basis_;
if (!highs_basis.useful) return;
const bool has_simplex_basis = ekk_instance_.status_.has_basis;
SimplexBasis& simplex_basis = ekk_instance_.basis_;
HighsLp& lp = model_.lp_;
assert(highs_basis.col_status.size() == static_cast<size_t>(lp.num_col_));
assert(highs_basis.row_status.size() == static_cast<size_t>(lp.num_row_));
HighsInt newNumCol = lp.num_col_ + ext_num_new_col;
HighsInt newNumTot = newNumCol + lp.num_row_;
highs_basis.col_status.resize(newNumCol);
if (has_simplex_basis) {
simplex_basis.nonbasicFlag_.resize(newNumTot);
simplex_basis.nonbasicMove_.resize(newNumTot);
for (HighsInt iRow = lp.num_row_ - 1; iRow >= 0; iRow--) {
HighsInt iCol = simplex_basis.basicIndex_[iRow];
if (iCol >= lp.num_col_) {
simplex_basis.basicIndex_[iRow] += ext_num_new_col;
}
simplex_basis.nonbasicFlag_[newNumCol + iRow] =
simplex_basis.nonbasicFlag_[lp.num_col_ + iRow];
simplex_basis.nonbasicMove_[newNumCol + iRow] =
simplex_basis.nonbasicMove_[lp.num_col_ + iRow];
}
}
for (HighsInt iCol = lp.num_col_; iCol < newNumCol; iCol++) {
double lower = lp.col_lower_[iCol];
double upper = lp.col_upper_[iCol];
HighsBasisStatus status = HighsBasisStatus::kNonbasic;
int8_t move = kIllegalMoveValue;
if (lower == upper) {
status = HighsBasisStatus::kLower;
move = kNonbasicMoveZe;
} else if (!highs_isInfinity(-lower)) {
if (!highs_isInfinity(upper)) {
if (fabs(lower) < fabs(upper)) {
status = HighsBasisStatus::kLower;
move = kNonbasicMoveUp;
} else {
status = HighsBasisStatus::kUpper;
move = kNonbasicMoveDn;
}
} else {
status = HighsBasisStatus::kLower;
move = kNonbasicMoveUp;
}
} else if (!highs_isInfinity(upper)) {
status = HighsBasisStatus::kUpper;
move = kNonbasicMoveDn;
} else {
status = HighsBasisStatus::kZero;
move = kNonbasicMoveZe;
}
assert(status != HighsBasisStatus::kNonbasic);
highs_basis.col_status[iCol] = status;
if (has_simplex_basis) {
assert(move != kIllegalMoveValue);
simplex_basis.nonbasicFlag_[iCol] = kNonbasicFlagTrue;
simplex_basis.nonbasicMove_[iCol] = move;
}
}
}
void Highs::appendBasicRowsToBasisInterface(const HighsInt ext_num_new_row) {
if (ext_num_new_row == 0) return;
HighsBasis& highs_basis = basis_;
if (!highs_basis.useful) return;
const bool has_simplex_basis = ekk_instance_.status_.has_basis;
SimplexBasis& simplex_basis = ekk_instance_.basis_;
HighsLp& lp = model_.lp_;
assert(highs_basis.col_status.size() == static_cast<size_t>(lp.num_col_));
assert(highs_basis.row_status.size() == static_cast<size_t>(lp.num_row_));
HighsInt newNumRow = lp.num_row_ + ext_num_new_row;
highs_basis.row_status.resize(newNumRow);
for (HighsInt iRow = lp.num_row_; iRow < newNumRow; iRow++)
highs_basis.row_status[iRow] = HighsBasisStatus::kBasic;
if (has_simplex_basis) {
HighsInt newNumTot = lp.num_col_ + newNumRow;
simplex_basis.nonbasicFlag_.resize(newNumTot);
simplex_basis.nonbasicMove_.resize(newNumTot);
simplex_basis.basicIndex_.resize(newNumRow);
for (HighsInt iRow = lp.num_row_; iRow < newNumRow; iRow++) {
simplex_basis.nonbasicFlag_[lp.num_col_ + iRow] = kNonbasicFlagFalse;
simplex_basis.nonbasicMove_[lp.num_col_ + iRow] = 0;
simplex_basis.basicIndex_[iRow] = lp.num_col_ + iRow;
}
}
}
HighsStatus Highs::getBasicVariablesInterface(HighsInt* basic_variables) {
HighsStatus return_status = HighsStatus::kOk;
HighsLp& lp = model_.lp_;
HighsInt num_row = lp.num_row_;
HighsInt num_col = lp.num_col_;
HighsSimplexStatus& ekk_status = ekk_instance_.status_;
if (num_row == 0) return return_status;
if (!basis_.valid) {
highsLogUser(options_.log_options, HighsLogType::kError,
"getBasicVariables called without a HiGHS basis\n");
return HighsStatus::kError;
}
if (!ekk_status.has_invert) {
HighsLpSolverObject solver_object(lp, basis_, solution_, info_,
ekk_instance_, callback_, options_,
timer_, sub_solver_call_time_);
const bool only_from_known_basis = true;
return_status = interpretCallStatus(
options_.log_options,
formSimplexLpBasisAndFactor(solver_object, only_from_known_basis),
return_status, "formSimplexLpBasisAndFactor");
if (return_status != HighsStatus::kOk) return return_status;
}
assert(ekk_status.has_invert);
for (HighsInt row = 0; row < num_row; row++) {
HighsInt var = ekk_instance_.basis_.basicIndex_[row];
if (var < num_col) {
basic_variables[row] = var;
} else {
basic_variables[row] = -(1 + var - num_col);
}
}
return return_status;
}
HighsStatus Highs::basisSolveInterface(const vector<double>& rhs,
double* solution_vector,
HighsInt* solution_num_nz,
HighsInt* solution_indices,
bool transpose) {
HighsStatus return_status = HighsStatus::kOk;
HighsLp& lp = model_.lp_;
HighsInt num_row = lp.num_row_;
if (num_row == 0) return return_status;
assert(ekk_instance_.status_.has_invert);
ekk_instance_.setNlaPointersForLpAndScale(lp);
assert(!lp.is_moved_);
HVector solve_vector;
solve_vector.setup(num_row);
solve_vector.clear();
HighsInt rhs_num_nz = 0;
for (HighsInt iRow = 0; iRow < num_row; iRow++) {
if (rhs[iRow]) {
solve_vector.index[rhs_num_nz++] = iRow;
solve_vector.array[iRow] = rhs[iRow];
}
}
solve_vector.count = rhs_num_nz;
const double expected_density = 1;
if (transpose) {
ekk_instance_.btran(solve_vector, expected_density);
} else {
ekk_instance_.ftran(solve_vector, expected_density);
}
if (solution_indices == NULL) {
if (solve_vector.count > num_row) {
for (HighsInt iRow = 0; iRow < num_row; iRow++) {
solution_vector[iRow] = solve_vector.array[iRow];
}
} else {
for (HighsInt iRow = 0; iRow < num_row; iRow++) solution_vector[iRow] = 0;
for (HighsInt iX = 0; iX < solve_vector.count; iX++) {
HighsInt iRow = solve_vector.index[iX];
solution_vector[iRow] = solve_vector.array[iRow];
}
}
} else {
if (solve_vector.count > num_row) {
solution_num_nz = 0;
for (HighsInt iRow = 0; iRow < num_row; iRow++) {
solution_vector[iRow] = 0;
if (solve_vector.array[iRow]) {
solution_vector[iRow] = solve_vector.array[iRow];
solution_indices[*solution_num_nz++] = iRow;
}
}
} else {
for (HighsInt iRow = 0; iRow < num_row; iRow++) solution_vector[iRow] = 0;
for (HighsInt iX = 0; iX < solve_vector.count; iX++) {
HighsInt iRow = solve_vector.index[iX];
solution_vector[iRow] = solve_vector.array[iRow];
solution_indices[iX] = iRow;
}
*solution_num_nz = solve_vector.count;
}
}
return HighsStatus::kOk;
}
void Highs::zeroIterationCounts() {
info_.simplex_iteration_count = 0;
info_.ipm_iteration_count = 0;
info_.crossover_iteration_count = 0;
info_.pdlp_iteration_count = 0;
info_.qp_iteration_count = 0;
}
HighsStatus Highs::getDualRayInterface(bool& has_dual_ray,
double* dual_ray_value) {
HighsStatus return_status = HighsStatus::kOk;
HighsLp& lp = model_.lp_;
HighsInt num_row = lp.num_row_;
if (num_row == 0) return return_status;
bool has_invert = ekk_instance_.status_.has_invert;
assert(!lp.is_moved_);
has_dual_ray = ekk_instance_.dual_ray_record_.index != kNoRayIndex;
std::vector<double> col_cost;
HighsHessian hessian;
bool solve_relaxation;
std::string presolve;
bool solve_feasibility_problem = false;
const bool is_qp = model_.isQp();
if (dual_ray_value != NULL) {
if (!has_dual_ray || !has_invert) {
if (this->model_status_ == HighsModelStatus::kOptimal) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Model status is optimal, so no dual ray is available\n");
return return_status;
}
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Solving LP to try to compute dual ray\n");
col_cost = lp.col_cost_;
if (is_qp) hessian = model_.hessian_;
this->getOptionValue("presolve", presolve);
this->getOptionValue("solve_relaxation", solve_relaxation);
solve_feasibility_problem = true;
std::vector<double> zero_costs;
zero_costs.assign(lp.num_col_, 0);
HighsRayRecord primal_ray_record =
this->ekk_instance_.primal_ray_record_.getRayRecord();
HighsStatus status =
this->changeColsCost(0, lp.num_col_ - 1, zero_costs.data());
assert(status == HighsStatus::kOk);
this->ekk_instance_.primal_ray_record_.setRayRecord(primal_ray_record);
if (is_qp) {
HighsHessian zero_hessian;
this->passHessian(zero_hessian);
}
this->setOptionValue("presolve", kHighsOffString);
this->setOptionValue("solve_relaxation", true);
HighsStatus call_status = this->run();
if (call_status != HighsStatus::kOk) return_status = call_status;
has_dual_ray = ekk_instance_.dual_ray_record_.index != kNoRayIndex;
has_invert = ekk_instance_.status_.has_invert;
assert(has_invert);
}
if (has_dual_ray) {
if (ekk_instance_.dual_ray_record_.value.size()) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Copying known dual ray\n");
for (HighsInt iRow = 0; iRow < num_row; iRow++)
dual_ray_value[iRow] = ekk_instance_.dual_ray_record_.value[iRow];
} else if (has_invert) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Solving linear system to compute dual ray\n");
vector<double> rhs;
HighsInt iRow = ekk_instance_.dual_ray_record_.index;
rhs.assign(num_row, 0);
rhs[iRow] = ekk_instance_.dual_ray_record_.sign;
HighsInt* dual_ray_num_nz = 0;
basisSolveInterface(rhs, dual_ray_value, dual_ray_num_nz, NULL, true);
ekk_instance_.dual_ray_record_.value.resize(num_row);
for (HighsInt iRow = 0; iRow < num_row; iRow++)
ekk_instance_.dual_ray_record_.value[iRow] = dual_ray_value[iRow];
} else {
assert(!has_invert);
highsLogUser(options_.log_options, HighsLogType::kError,
"No LP invertible representation to compute dual ray\n");
return_status = HighsStatus::kError;
}
} else {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"No dual ray found\n");
return_status = HighsStatus::kOk;
}
}
if (solve_feasibility_problem) {
lp.col_cost_ = col_cost;
if (is_qp) model_.hessian_ = hessian;
this->setOptionValue("presolve", presolve);
this->setOptionValue("solve_relaxation", solve_relaxation);
this->info_.primal_solution_status = SolutionStatus::kSolutionStatusNone;
this->info_.dual_solution_status = SolutionStatus::kSolutionStatusNone;
this->info_.objective_function_value = 0;
this->info_.invalidateDualKkt();
if (has_dual_ray) {
assert(this->info_.num_primal_infeasibilities > 0);
assert(this->model_status_ == HighsModelStatus::kInfeasible);
} else {
this->info_.invalidatePrimalKkt();
this->model_status_ = HighsModelStatus::kNotset;
}
}
return return_status;
}
HighsStatus Highs::getPrimalRayInterface(bool& has_primal_ray,
double* primal_ray_value) {
HighsStatus return_status = HighsStatus::kOk;
HighsLp& lp = model_.lp_;
HighsInt num_row = lp.num_row_;
HighsInt num_col = lp.num_col_;
if (num_row == 0) return return_status;
if (model_.isQp()) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Cannot find primal ray for unbounded QP\n");
return HighsStatus::kError;
}
bool has_invert = ekk_instance_.status_.has_invert;
assert(!lp.is_moved_);
has_primal_ray = ekk_instance_.primal_ray_record_.index != kNoRayIndex;
std::string presolve;
bool solve_relaxation;
bool allow_unbounded_or_infeasible;
bool solve_unboundedness_problem = false;
if (primal_ray_value != NULL) {
if (!has_primal_ray || !has_invert) {
if (this->model_status_ == HighsModelStatus::kOptimal) {
highsLogUser(
options_.log_options, HighsLogType::kInfo,
"Model status is optimal, so no primal ray is available\n");
return return_status;
}
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Solving LP to try to compute primal ray\n");
this->getOptionValue("presolve", presolve);
this->getOptionValue("solve_relaxation", solve_relaxation);
this->getOptionValue("allow_unbounded_or_infeasible",
allow_unbounded_or_infeasible);
solve_unboundedness_problem = true;
this->setOptionValue("presolve", kHighsOffString);
this->setOptionValue("solve_relaxation", true);
this->setOptionValue("allow_unbounded_or_infeasible", false);
HighsStatus call_status = this->run();
if (call_status != HighsStatus::kOk) return_status = call_status;
has_primal_ray = ekk_instance_.primal_ray_record_.index != kNoRayIndex;
has_invert = ekk_instance_.status_.has_invert;
assert(has_invert);
}
if (has_primal_ray) {
if (ekk_instance_.primal_ray_record_.value.size()) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Copying known primal ray\n");
for (HighsInt iCol = 0; iCol < num_col; iCol++)
primal_ray_value[iCol] = ekk_instance_.primal_ray_record_.value[iCol];
return return_status;
} else if (has_invert) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Solving linear system to compute primal ray\n");
HighsInt col = ekk_instance_.primal_ray_record_.index;
assert(ekk_instance_.basis_.nonbasicFlag_[col] == kNonbasicFlagTrue);
vector<double> rhs;
vector<double> column;
column.assign(num_row, 0);
rhs.assign(num_row, 0);
lp.ensureColwise();
HighsInt primal_ray_sign = ekk_instance_.primal_ray_record_.sign;
if (col < num_col) {
for (HighsInt iEl = lp.a_matrix_.start_[col];
iEl < lp.a_matrix_.start_[col + 1]; iEl++)
rhs[lp.a_matrix_.index_[iEl]] =
primal_ray_sign * lp.a_matrix_.value_[iEl];
} else {
rhs[col - num_col] = primal_ray_sign;
}
HighsInt* column_num_nz = 0;
basisSolveInterface(rhs, column.data(), column_num_nz, NULL, false);
for (HighsInt iCol = 0; iCol < num_col; iCol++)
primal_ray_value[iCol] = 0;
for (HighsInt iRow = 0; iRow < num_row; iRow++) {
HighsInt iCol = ekk_instance_.basis_.basicIndex_[iRow];
if (iCol < num_col) primal_ray_value[iCol] = column[iRow];
}
if (col < num_col) primal_ray_value[col] = -primal_ray_sign;
ekk_instance_.primal_ray_record_.value.resize(num_col);
for (HighsInt iCol = 0; iCol < num_col; iCol++)
ekk_instance_.primal_ray_record_.value[iCol] = primal_ray_value[iCol];
}
} else {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"No primal ray found\n");
return_status = HighsStatus::kOk;
}
}
const bool is_mip = this->model_.isMip();
if (solve_unboundedness_problem) {
if (is_mip) {
this->info_.dual_solution_status = SolutionStatus::kSolutionStatusNone;
this->info_.invalidateDualKkt();
}
this->setOptionValue("presolve", presolve);
this->setOptionValue("solve_relaxation", solve_relaxation);
this->setOptionValue("allow_unbounded_or_infeasible",
allow_unbounded_or_infeasible);
if (has_primal_ray) {
assert(is_mip || this->info_.num_dual_infeasibilities > 0);
assert(this->model_status_ == HighsModelStatus::kUnbounded);
}
}
return return_status;
}
HighsStatus Highs::getRangingInterface() {
HighsLpSolverObject solver_object(model_.lp_, basis_, solution_, info_,
ekk_instance_, callback_, options_, timer_,
sub_solver_call_time_);
solver_object.model_status_ = model_status_;
return getRangingData(this->ranging_, solver_object);
}
HighsStatus Highs::getIisInterfaceReturn(const HighsStatus return_status) {
if (return_status == HighsStatus::kError) return return_status;
assert(this->iis_.valid_);
if (this->iis_.status_ >= kIisModelStatusInfeasible)
this->model_status_ = HighsModelStatus::kInfeasible;
HighsLp& lp = this->model_.lp_;
HighsLp& iis_lp = this->iis_.model_.lp_;
const bool has_is =
this->iis_.col_index_.size() || this->iis_.row_index_.size();
const bool has_iis = this->iis_.status_ == kIisModelStatusIrreducible;
if (has_iis) assert(has_is);
if (has_is) {
this->iis_.setLp(lp);
bool lp_data_ok = this->iis_.lpDataOk(lp, this->options_);
assert(lp_data_ok);
if (!lp_data_ok) return HighsStatus::kError;
bool lp_ok = this->iis_.lpOk(this->options_);
if (!lp_ok) return HighsStatus::kError;
} else {
assert(this->iis_.status_ <= kIisModelStatusReducible);
assert(!this->iis_.col_bound_.size());
assert(!this->iis_.row_bound_.size());
}
this->iis_.setStatus(lp);
bool index_status_ok = this->iis_.indexStatusOk(lp);
assert(index_status_ok);
if (!index_status_ok) return HighsStatus::kError;
return return_status;
}
HighsStatus Highs::getIisInterface() {
const HighsLp& lp = model_.lp_;
if (this->model_status_ == HighsModelStatus::kOptimal ||
this->model_status_ == HighsModelStatus::kUnbounded) {
highsLogUser(
options_.log_options, HighsLogType::kInfo,
"Calling Highs::getIis for a model that is known to be feasible\n");
this->iis_.clear();
this->iis_.valid_ = true;
this->iis_.status_ = kIisModelStatusFeasible;
return this->getIisInterfaceReturn(HighsStatus::kOk);
}
HighsStatus return_status = HighsStatus::kOk;
if (this->model_status_ != HighsModelStatus::kNotset &&
this->model_status_ != HighsModelStatus::kInfeasible) {
return_status = HighsStatus::kWarning;
highsLogUser(options_.log_options, HighsLogType::kWarning,
"Calling Highs::getIis for a model with status %s\n",
this->modelStatusToString(this->model_status_).c_str());
}
if (this->iis_.valid_) return this->getIisInterfaceReturn(HighsStatus::kOk);
this->iis_.clear();
if (this->iis_.trivial(lp, options_))
return this->getIisInterfaceReturn(HighsStatus::kOk);
HighsInt num_row = lp.num_row_;
if (num_row == 0) {
this->iis_.valid_ = true;
return this->getIisInterfaceReturn(HighsStatus::kOk);
}
if (this->iis_.rowValueBounds(lp, options_))
return this->getIisInterfaceReturn(HighsStatus::kOk);
if (options_.iis_strategy == kIisStrategyLight)
return this->getIisInterfaceReturn(HighsStatus::kOk);
bool ray_option =
false;
const bool lp_option = options_.iis_strategy >= kIisStrategyFromLp;
if (this->model_status_ == HighsModelStatus::kInfeasible && ray_option &&
!ekk_instance_.status_.has_invert) {
std::string presolve = options_.presolve;
options_.presolve = kHighsOffString;
HighsIisInfo iis_info;
iis_info.simplex_time = -this->getRunTime();
iis_info.simplex_iterations = -info_.simplex_iteration_count;
HighsStatus run_status = this->optimizeModel();
options_.presolve = presolve;
if (run_status != HighsStatus::kOk) return run_status;
iis_info.simplex_time += this->getRunTime();
iis_info.simplex_iterations += -info_.simplex_iteration_count;
this->iis_.info_.push_back(iis_info);
if (this->model_status_ != HighsModelStatus::kInfeasible) {
highsLogUser(
options_.log_options, HighsLogType::kError,
"Model status has switched from %s to %s when solving without "
"presolve\n",
this->modelStatusToString(HighsModelStatus::kInfeasible).c_str(),
this->modelStatusToString(this->model_status_).c_str());
return HighsStatus::kError;
}
}
const bool has_dual_ray = ekk_instance_.dual_ray_record_.index != kNoRayIndex;
if (ray_option && !has_dual_ray)
highsLogUser(options_.log_options, HighsLogType::kWarning,
"No known dual ray from which to compute IIS\n");
ray_option = ray_option && has_dual_ray;
if (ray_option) {
assert(has_dual_ray);
assert(ekk_instance_.status_.has_invert);
assert(!lp.is_moved_);
std::vector<double> rhs;
HighsInt iRow = ekk_instance_.dual_ray_record_.index;
rhs.assign(num_row, 0);
rhs[iRow] = 1;
std::vector<double> dual_ray_value(num_row);
HighsInt* dual_ray_num_nz = 0;
basisSolveInterface(rhs, dual_ray_value.data(), dual_ray_num_nz, NULL,
true);
for (HighsInt iRow = 0; iRow < lp.num_row_; iRow++)
if (dual_ray_value[iRow]) this->iis_.row_index_.push_back(iRow);
} else if (lp_option) {
this->invalidateSolverData();
HighsLp check_lp_before = this->model_.lp_;
HighsStatus return_status = this->elasticityFilter(-1.0, -1.0, 1.0, nullptr,
nullptr, nullptr, true);
HighsLp check_lp_after = this->model_.lp_;
assert(check_lp_before.equalVectors(check_lp_after));
assert(check_lp_before.a_matrix_.equivalent(check_lp_after.a_matrix_));
if (return_status != HighsStatus::kOk) return return_status;
}
if (!ray_option && !lp_option)
return this->getIisInterfaceReturn(HighsStatus::kOk);
HighsIis iis = this->iis_;
this->invalidateSolverData();
this->iis_ = iis;
return_status = HighsStatus::kOk;
if (this->iis_.row_index_.size() == 0) {
assert(this->iis_.valid_);
assert(this->iis_.status_ == kIisModelStatusFeasible);
return this->getIisInterfaceReturn(return_status);
}
assert(this->iis_.row_index_.size());
this->iis_.status_ = kIisModelStatusReducible;
if (!(kIisStrategyIrreducible & this->options_.iis_strategy))
return this->getIisInterfaceReturn(return_status);
model_.lp_.a_matrix_.ensureColwise();
return_status = this->iis_.deduce(lp, options_, basis_);
const HighsInt num_lp_solved = this->iis_.info_.size();
double min_time = kHighsInf;
double sum_time = 0;
double max_time = 0;
HighsInt min_iterations = kHighsIInf;
HighsInt sum_iterations = 0;
HighsInt max_iterations = 0;
for (HighsInt iX = 0; iX < num_lp_solved; iX++) {
double time = this->iis_.info_[iX].simplex_time;
HighsInt iterations = this->iis_.info_[iX].simplex_iterations;
min_time = std::min(time, min_time);
sum_time += time;
max_time = std::max(time, max_time);
min_iterations = std::min(iterations, min_iterations);
sum_iterations += iterations;
max_iterations = std::max(iterations, max_iterations);
}
highsLogUser(options_.log_options, HighsLogType::kInfo,
" %d cols, %d rows, %d LPs solved"
" (min / average / max) iteration count (%6d / %6.2g / % 6d)"
" and time (%6.2f / %6.2f / % 6.2f) \n",
int(this->iis_.col_index_.size()),
int(this->iis_.row_index_.size()), int(num_lp_solved),
int(min_iterations),
num_lp_solved > 0 ? (1.0 * sum_iterations) / num_lp_solved : 0,
int(max_iterations), min_time,
num_lp_solved > 0 ? sum_time / num_lp_solved : 0, max_time);
return this->getIisInterfaceReturn(return_status);
}
HighsStatus Highs::elasticityFilterReturn(
const HighsStatus return_status, const std::string& original_model_name,
const HighsInt original_num_col, const HighsInt original_num_row,
const std::vector<double>& original_col_cost,
const std::vector<double>& original_col_lower,
const std::vector<double>& original_col_upper,
const std::vector<HighsVarType>& original_integrality) {
const HighsLp& lp = this->model_.lp_;
const bool iis_valid = this->iis_.valid_;
HighsIis iis;
if (iis_valid) {
iis = this->iis_;
} else {
assert(return_status != HighsStatus::kOk);
}
double objective_function_value = info_.objective_function_value;
HighsStatus run_status;
run_status = this->deleteRows(original_num_row, lp.num_row_ - 1);
assert(run_status == HighsStatus::kOk);
run_status = this->deleteCols(original_num_col, lp.num_col_ - 1);
assert(run_status == HighsStatus::kOk);
basis_.valid = false;
run_status =
this->changeColsCost(0, original_num_col - 1, original_col_cost.data());
assert(run_status == HighsStatus::kOk);
run_status =
this->changeColsBounds(0, original_num_col - 1, original_col_lower.data(),
original_col_upper.data());
assert(run_status == HighsStatus::kOk);
if (kIisStrategyRelaxation & this->options_.iis_strategy &&
original_integrality.size()) {
this->changeColsIntegrality(0, original_num_col - 1,
original_integrality.data());
assert(run_status == HighsStatus::kOk);
}
this->passModelName(original_model_name);
assert(lp.num_col_ == original_num_col);
assert(lp.num_row_ == original_num_row);
if (return_status == HighsStatus::kOk) {
this->model_.lp_.a_matrix_.productQuad(this->solution_.row_value,
this->solution_.col_value);
this->solution_.value_valid = true;
info_.objective_function_value = objective_function_value;
getKktFailures(options_, model_, solution_, basis_, info_);
info_.valid = true;
}
if (iis_valid) {
if (iis.status_ == kIisModelStatusFeasible)
this->model_status_ = HighsModelStatus::kNotset;
this->iis_ = iis;
} else {
this->model_status_ = HighsModelStatus::kNotset;
}
return return_status;
}
HighsStatus Highs::elasticityFilter(const double global_lower_penalty,
const double global_upper_penalty,
const double global_rhs_penalty,
const double* local_lower_penalty,
const double* local_upper_penalty,
const double* local_rhs_penalty,
const bool get_iis) {
std::vector<HighsInt> col_of_ecol;
std::vector<HighsInt> row_of_ecol;
std::vector<bool> bound_of_row_of_ecol_is_lower;
std::vector<bool> bound_of_col_of_ecol_is_lower;
std::vector<double> erow_lower;
std::vector<double> erow_upper;
std::vector<HighsInt> erow_start;
std::vector<HighsInt> erow_index;
std::vector<double> erow_value;
std::vector<std::string> ecol_name;
std::vector<std::string> erow_name;
std::vector<double> ecol_cost;
std::vector<double> ecol_lower;
std::vector<double> ecol_upper;
const HighsLp& lp = this->model_.lp_;
HighsInt evar_ix = lp.num_col_;
HighsStatus run_status;
const bool write_model = false;
const std::string original_model_name = lp.model_name_;
const HighsInt original_num_col = lp.num_col_;
const HighsInt original_num_row = lp.num_row_;
const std::vector<double> original_col_cost = lp.col_cost_;
const std::vector<double> original_col_lower = lp.col_lower_;
const std::vector<double> original_col_upper = lp.col_upper_;
const std::vector<HighsVarType> original_integrality = lp.integrality_;
this->passModelName(original_model_name + "_elastic");
std::vector<double> zero_costs;
zero_costs.assign(original_num_col, 0);
run_status = this->changeColsCost(0, lp.num_col_ - 1, zero_costs.data());
assert(run_status == HighsStatus::kOk);
const bool mip_relaxation_iis =
get_iis && lp.isMip() &&
kIisStrategyRelaxation & this->options_.iis_strategy;
if (mip_relaxation_iis) {
std::vector<HighsVarType> all_continuous;
all_continuous.assign(original_num_col, HighsVarType::kContinuous);
run_status =
this->changeColsIntegrality(0, lp.num_col_ - 1, all_continuous.data());
assert(run_status == HighsStatus::kOk);
}
const bool has_local_lower_penalty = local_lower_penalty;
const bool has_global_elastic_lower = global_lower_penalty >= 0;
const bool has_elastic_lower =
has_local_lower_penalty || has_global_elastic_lower;
const bool has_local_upper_penalty = local_upper_penalty;
const bool has_global_elastic_upper = global_upper_penalty >= 0;
const bool has_elastic_upper =
has_local_upper_penalty || has_global_elastic_upper;
const bool has_elastic_columns = has_elastic_lower || has_elastic_upper;
const bool has_local_rhs_penalty = local_rhs_penalty;
const bool has_global_elastic_rhs = global_rhs_penalty >= 0;
const bool has_elastic_rows = has_local_rhs_penalty || has_global_elastic_rhs;
assert(has_elastic_columns || has_elastic_rows);
const HighsInt col_ecol_offset = lp.num_col_;
if (has_elastic_columns) {
std::vector<double> col_lower;
std::vector<double> col_upper;
const bool has_col_names = lp.col_names_.size() > 0;
erow_start.push_back(0);
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
const double lower = lp.col_lower_[iCol];
const double upper = lp.col_upper_[iCol];
col_lower.push_back(lower);
col_upper.push_back(upper);
if (lower <= -kHighsInf && upper >= kHighsInf) continue;
const double lower_penalty = has_local_lower_penalty
? local_lower_penalty[iCol]
: global_lower_penalty;
if (lower_penalty < 0 && upper >= kHighsInf) continue;
const double upper_penalty = has_local_upper_penalty
? local_upper_penalty[iCol]
: global_upper_penalty;
if (lower <= -kHighsInf && upper_penalty < 0) continue;
erow_lower.push_back(lower);
erow_upper.push_back(upper);
if (has_col_names)
erow_name.push_back("row_" + std::to_string(iCol) + "_" +
lp.col_names_[iCol] + "_erow");
erow_index.push_back(iCol);
erow_value.push_back(1);
if (lower > -kHighsInf && lower_penalty >= 0) {
col_of_ecol.push_back(iCol);
if (has_col_names)
ecol_name.push_back("col_" + std::to_string(iCol) + "_" +
lp.col_names_[iCol] + "_lower");
bound_of_col_of_ecol_is_lower.push_back(true);
col_lower[iCol] = -kHighsInf;
erow_index.push_back(evar_ix);
erow_value.push_back(1);
ecol_cost.push_back(lower_penalty);
evar_ix++;
}
if (upper < kHighsInf && upper_penalty >= 0) {
col_of_ecol.push_back(iCol);
if (has_col_names)
ecol_name.push_back("col_" + std::to_string(iCol) + "_" +
lp.col_names_[iCol] + "_upper");
bound_of_col_of_ecol_is_lower.push_back(false);
col_upper[iCol] = kHighsInf;
erow_index.push_back(evar_ix);
erow_value.push_back(-1);
ecol_cost.push_back(upper_penalty);
evar_ix++;
}
erow_start.push_back(erow_index.size());
HighsInt row_nz =
erow_start[erow_start.size() - 1] - erow_start[erow_start.size() - 2];
assert(row_nz == 2 || row_nz == 3);
}
HighsInt num_new_col = col_of_ecol.size();
HighsInt num_new_row = erow_start.size() - 1;
HighsInt num_new_nz = erow_start[num_new_row];
if (kIisDevReport)
printf(
"Elasticity filter: For columns there are %d variables and %d "
"constraints\n",
int(num_new_col), int(num_new_row));
assert(col_lower.size() == static_cast<size_t>(lp.num_col_));
assert(col_upper.size() == static_cast<size_t>(lp.num_col_));
run_status = this->changeColsBounds(0, lp.num_col_ - 1, col_lower.data(),
col_upper.data());
assert(run_status == HighsStatus::kOk);
ecol_lower.assign(num_new_col, 0);
ecol_upper.assign(num_new_col, kHighsInf);
run_status = this->addCols(num_new_col, ecol_cost.data(), ecol_lower.data(),
ecol_upper.data(), 0, nullptr, nullptr, nullptr);
assert(run_status == HighsStatus::kOk);
assert(erow_start.size() == static_cast<size_t>(num_new_row + 1));
assert(erow_index.size() == static_cast<size_t>(num_new_nz));
assert(erow_value.size() == static_cast<size_t>(num_new_nz));
run_status = this->addRows(num_new_row, erow_lower.data(),
erow_upper.data(), num_new_nz, erow_start.data(),
erow_index.data(), erow_value.data());
assert(run_status == HighsStatus::kOk);
if (has_col_names) {
for (HighsInt iCol = 0; iCol < num_new_col; iCol++)
this->passColName(col_ecol_offset + iCol, ecol_name[iCol]);
for (HighsInt iRow = 0; iRow < num_new_row; iRow++)
this->passRowName(original_num_row + iRow, erow_name[iRow]);
}
assert(ecol_cost.size() == static_cast<size_t>(num_new_col));
assert(ecol_lower.size() == static_cast<size_t>(num_new_col));
assert(ecol_upper.size() == static_cast<size_t>(num_new_col));
if (write_model) {
printf("\nAfter adding %d e-rows\n=============\n", int(num_new_col));
bool output_flag;
run_status = this->getOptionValue("output_flag", output_flag);
this->setOptionValue("output_flag", true);
this->writeModel("");
this->setOptionValue("output_flag", output_flag);
}
}
const HighsInt row_ecol_offset = lp.num_col_;
if (has_elastic_rows) {
ecol_name.clear();
ecol_cost.clear();
std::vector<HighsInt> ecol_start;
std::vector<HighsInt> ecol_index;
std::vector<double> ecol_value;
ecol_start.push_back(0);
const bool has_row_names = lp.row_names_.size() > 0;
for (HighsInt iRow = 0; iRow < original_num_row; iRow++) {
const double penalty =
has_local_rhs_penalty ? local_rhs_penalty[iRow] : global_rhs_penalty;
if (penalty < 0) continue;
const double lower = lp.row_lower_[iRow];
const double upper = lp.row_upper_[iRow];
if (lower > -kHighsInf) {
row_of_ecol.push_back(iRow);
if (has_row_names)
ecol_name.push_back("row_" + std::to_string(iRow) + "_" +
lp.row_names_[iRow] + "_lower");
bound_of_row_of_ecol_is_lower.push_back(true);
ecol_index.push_back(iRow);
ecol_value.push_back(1);
ecol_start.push_back(ecol_index.size());
ecol_cost.push_back(penalty);
evar_ix++;
}
if (upper < kHighsInf) {
row_of_ecol.push_back(iRow);
if (has_row_names)
ecol_name.push_back("row_" + std::to_string(iRow) + "_" +
lp.row_names_[iRow] + "_upper");
bound_of_row_of_ecol_is_lower.push_back(false);
ecol_index.push_back(iRow);
ecol_value.push_back(-1);
ecol_start.push_back(ecol_index.size());
ecol_cost.push_back(penalty);
evar_ix++;
}
}
HighsInt num_new_col = ecol_start.size() - 1;
HighsInt num_new_nz = ecol_start[num_new_col];
ecol_lower.assign(num_new_col, 0);
ecol_upper.assign(num_new_col, kHighsInf);
assert(ecol_cost.size() == static_cast<size_t>(num_new_col));
assert(ecol_lower.size() == static_cast<size_t>(num_new_col));
assert(ecol_upper.size() == static_cast<size_t>(num_new_col));
assert(ecol_start.size() == static_cast<size_t>(num_new_col + 1));
assert(ecol_index.size() == static_cast<size_t>(num_new_nz));
assert(ecol_value.size() == static_cast<size_t>(num_new_nz));
run_status = this->addCols(num_new_col, ecol_cost.data(), ecol_lower.data(),
ecol_upper.data(), num_new_nz, ecol_start.data(),
ecol_index.data(), ecol_value.data());
assert(run_status == HighsStatus::kOk);
if (has_row_names) {
for (HighsInt iCol = 0; iCol < num_new_col; iCol++)
this->passColName(row_ecol_offset + iCol, ecol_name[iCol]);
}
if (write_model) {
bool output_flag;
printf("\nAfter adding %d e-cols\n=============\n", int(num_new_col));
run_status = this->getOptionValue("output_flag", output_flag);
this->setOptionValue("output_flag", true);
this->writeModel("");
this->setOptionValue("output_flag", output_flag);
}
}
if (write_model) this->writeModel("elastic.mps");
HighsIis& iis = this->iis_;
iis.clear();
auto solveLp = [&]() -> HighsStatus {
HighsIisInfo iis_info;
iis_info.simplex_time = -this->getRunTime();
iis_info.simplex_iterations = -info_.simplex_iteration_count;
run_status = this->optimizeModel();
assert(run_status == HighsStatus::kOk);
if (run_status != HighsStatus::kOk) return run_status;
iis_info.simplex_time += this->getRunTime();
iis_info.simplex_iterations += info_.simplex_iteration_count;
iis.info_.push_back(iis_info);
return run_status;
};
run_status = solveLp();
this->writeSolution("", 1);
if (run_status != HighsStatus::kOk)
return elasticityFilterReturn(run_status, original_model_name,
original_num_col, original_num_row,
original_col_cost, original_col_lower,
original_col_upper, original_integrality);
if (kIisDevReport) this->writeSolution("", kSolutionStylePretty);
assert(this->model_status_ == HighsModelStatus::kOptimal ||
this->model_status_ == HighsModelStatus::kUnbounded);
this->iis_.valid_ = true;
this->iis_.status_ = this->info_.objective_function_value > 0
? kIisModelStatusInfeasible
: kIisModelStatusFeasible;
if (!get_iis)
return elasticityFilterReturn(HighsStatus::kOk, original_model_name,
original_num_col, original_num_row,
original_col_cost, original_col_lower,
original_col_upper, original_integrality);
assert(!has_elastic_columns);
assert(has_elastic_rows);
assert(original_num_row == lp.num_row_);
const HighsSolution& solution = this->getSolution();
HighsInt loop_k = 0;
for (;;) {
if (kIisDevReport)
printf("\nElasticity filter pass %d\n==============\n", int(loop_k));
HighsInt num_fixed = 0;
if (has_elastic_columns) {
for (size_t eCol = 0; eCol < col_of_ecol.size(); eCol++) {
HighsInt iCol = col_of_ecol[eCol];
if (solution.col_value[col_ecol_offset + eCol] >
this->options_.primal_feasibility_tolerance) {
if (kIisDevReport)
printf(
"E-col %2d (column %2d) corresponds to column %2d with %s "
"bound "
"%11.4g "
"and has solution value %11.4g\n",
int(eCol), int(col_ecol_offset + eCol), int(iCol),
bound_of_col_of_ecol_is_lower[eCol] ? "lower" : "upper",
bound_of_col_of_ecol_is_lower[eCol] ? lp.col_lower_[iCol]
: lp.col_upper_[iCol],
solution.col_value[col_ecol_offset + eCol]);
this->changeColBounds(col_ecol_offset + eCol, 0, 0);
num_fixed++;
}
}
}
if (has_elastic_rows) {
for (size_t eCol = 0; eCol < row_of_ecol.size(); eCol++) {
HighsInt iRow = row_of_ecol[eCol];
if (solution.col_value[row_ecol_offset + eCol] >
this->options_.primal_feasibility_tolerance) {
if (kIisDevReport)
printf(
"E-row %2d (column %2d) corresponds to row %2d with %s "
"bound "
"%11.4g "
"and has solution value %11.4g\n",
int(eCol), int(row_ecol_offset + eCol), int(iRow),
bound_of_row_of_ecol_is_lower[eCol] ? "lower" : "upper",
bound_of_row_of_ecol_is_lower[eCol] ? lp.row_lower_[iRow]
: lp.row_upper_[iRow],
solution.col_value[row_ecol_offset + eCol]);
this->changeColBounds(row_ecol_offset + eCol, 0, 0);
num_fixed++;
}
}
}
if (num_fixed == 0) {
this->iis_.status_ = kIisModelStatusFeasible;
break;
}
HighsStatus run_status = solveLp();
if (run_status != HighsStatus::kOk)
return elasticityFilterReturn(run_status, original_model_name,
original_num_col, original_num_row,
original_col_cost, original_col_lower,
original_col_upper, original_integrality);
if (kIisDevReport) this->writeSolution("", kSolutionStylePretty);
HighsModelStatus model_status = this->getModelStatus();
if (model_status == HighsModelStatus::kInfeasible) break;
loop_k++;
}
HighsInt num_enforced_col_ecol = 0;
HighsInt num_enforced_row_ecol = 0;
if (has_elastic_columns) {
for (size_t eCol = 0; eCol < col_of_ecol.size(); eCol++) {
HighsInt iCol = col_of_ecol[eCol];
if (lp.col_upper_[col_ecol_offset + eCol] == 0) {
num_enforced_col_ecol++;
if (kIisDevReport)
printf(
"Col e-col %2d (column %2d) corresponds to column %2d with %s "
"bound "
"%11.4g "
"and is enforced\n",
int(eCol), int(col_ecol_offset + eCol), int(iCol),
bound_of_col_of_ecol_is_lower[eCol] ? "lower" : "upper",
bound_of_col_of_ecol_is_lower[eCol] ? lp.col_lower_[iCol]
: lp.col_upper_[iCol]);
}
}
}
if (has_elastic_rows) {
for (size_t eCol = 0; eCol < row_of_ecol.size(); eCol++) {
HighsInt iRow = row_of_ecol[eCol];
if (lp.col_upper_[row_ecol_offset + eCol] == 0) {
num_enforced_row_ecol++;
iis.row_index_.push_back(iRow);
if (kIisDevReport)
printf(
"Row e-col %2d (column %2d) corresponds to row %2d with %s "
"bound "
"%11.4g and is enforced\n",
int(eCol), int(row_ecol_offset + eCol), int(iRow),
bound_of_row_of_ecol_is_lower[eCol] ? "lower" : "upper",
bound_of_row_of_ecol_is_lower[eCol] ? lp.row_lower_[iRow]
: lp.row_upper_[iRow]);
}
}
}
HighsInt num_iis_row = iis.row_index_.size();
if (iis.status_ == kIisModelStatusFeasible) {
assert(num_enforced_col_ecol == 0 && num_enforced_row_ecol == 0);
assert(num_iis_row == 0);
}
highsLogUser(
options_.log_options, HighsLogType::kInfo,
"Elasticity filter after %d passes enforces bounds on %d cols and %d "
"rows\n",
int(loop_k), int(num_enforced_col_ecol), int(num_enforced_row_ecol));
if (kIisDevReport)
printf(
"\nElasticity filter after %d passes enforces bounds on %d cols and %d "
"rows\n",
int(loop_k), int(num_enforced_col_ecol), int(num_enforced_row_ecol));
iis.valid_ = true;
iis.strategy_ = this->options_.iis_strategy;
if (iis.status_ == kIisModelStatusFeasible)
return elasticityFilterReturn(HighsStatus::kOk, original_model_name,
original_num_col, original_num_row,
original_col_cost, original_col_lower,
original_col_upper, original_integrality);
assert(this->iis_.status_ > kIisModelStatusFeasible);
assert(num_iis_row > 0);
iis.status_ = kIisModelStatusReducible;
std::vector<HighsInt> in_row_index(original_num_row, -1);
for (HighsInt iX = 0; iX < num_iis_row; iX++)
in_row_index[iis.row_index_[iX]] = iX;
std::vector<bool> nonzero_in_row_index(original_num_col, false);
if (lp.a_matrix_.isColwise()) {
for (HighsInt iCol = 0; iCol < original_num_col; iCol++) {
for (HighsInt iEl = lp.a_matrix_.start_[iCol];
iEl < lp.a_matrix_.start_[iCol + 1]; iEl++)
if (in_row_index[lp.a_matrix_.index_[iEl]] >= 0)
nonzero_in_row_index[iCol] = true;
}
} else {
for (HighsInt iX = 0; iX < num_iis_row; iX++) {
HighsInt iRow = iis.row_index_[iX];
for (HighsInt iEl = lp.a_matrix_.start_[iRow];
iEl < lp.a_matrix_.start_[iRow + 1]; iEl++)
nonzero_in_row_index[lp.a_matrix_.index_[iEl]] = true;
}
}
assert(iis.col_index_.size() == 0);
for (HighsInt iCol = 0; iCol < original_num_col; iCol++) {
if (nonzero_in_row_index[iCol]) {
iis.col_index_.push_back(iCol);
if (lp.col_lower_[iCol] > -kHighsInf) {
if (lp.col_upper_[iCol] < kHighsInf) {
iis.col_bound_.push_back(kIisBoundStatusBoxed);
} else {
iis.col_bound_.push_back(kIisBoundStatusLower);
}
} else {
if (lp.col_upper_[iCol] < kHighsInf) {
iis.col_bound_.push_back(kIisBoundStatusUpper);
} else {
iis.col_bound_.push_back(kIisBoundStatusFree);
}
}
}
}
iis.row_bound_.assign(num_iis_row, -1);
for (size_t eCol = 0; eCol < row_of_ecol.size(); eCol++) {
HighsInt iRow = row_of_ecol[eCol];
if (lp.col_upper_[row_ecol_offset + eCol] == 0) {
HighsInt iX = in_row_index[iRow];
assert(iis.row_bound_[iX] == -1);
iis.row_bound_[iX] = bound_of_row_of_ecol_is_lower[eCol]
? kIisBoundStatusLower
: kIisBoundStatusUpper;
}
}
assert(iis.row_bound_.size() == iis.row_index_.size());
for (HighsInt iX = 0; iX < num_iis_row; iX++) assert(iis.row_bound_[iX] >= 0);
return elasticityFilterReturn(HighsStatus::kOk, original_model_name,
original_num_col, original_num_row,
original_col_cost, original_col_lower,
original_col_upper, original_integrality);
}
HighsStatus Highs::extractIis(HighsInt& num_iis_col, HighsInt& num_iis_row,
HighsInt* iis_col_index, HighsInt* iis_row_index,
HighsInt* iis_col_bound,
HighsInt* iis_row_bound) {
assert(this->iis_.valid_);
num_iis_col = this->iis_.col_index_.size();
num_iis_row = this->iis_.row_index_.size();
if (iis_col_index || iis_col_bound) {
for (HighsInt iCol = 0; iCol < num_iis_col; iCol++) {
if (iis_col_index) iis_col_index[iCol] = this->iis_.col_index_[iCol];
if (iis_col_bound) iis_col_bound[iCol] = this->iis_.col_bound_[iCol];
}
}
if (iis_row_index || iis_row_bound) {
for (HighsInt iRow = 0; iRow < num_iis_row; iRow++) {
if (iis_row_index) iis_row_index[iRow] = this->iis_.row_index_[iRow];
if (iis_row_bound) iis_row_bound[iRow] = this->iis_.row_bound_[iRow];
}
}
return HighsStatus::kOk;
}
bool Highs::aFormatOk(const HighsInt num_nz, const HighsInt format) {
if (!num_nz) return true;
const bool ok_format = format == (HighsInt)MatrixFormat::kColwise ||
format == (HighsInt)MatrixFormat::kRowwise;
assert(ok_format);
if (!ok_format)
highsLogUser(
options_.log_options, HighsLogType::kError,
"Non-empty Constraint matrix has illegal format = %" HIGHSINT_FORMAT
"\n",
format);
return ok_format;
}
bool Highs::qFormatOk(const HighsInt num_nz, const HighsInt format) {
if (!num_nz) return true;
const bool ok_format = format == (HighsInt)HessianFormat::kTriangular;
assert(ok_format);
if (!ok_format)
highsLogUser(
options_.log_options, HighsLogType::kError,
"Non-empty Hessian matrix has illegal format = %" HIGHSINT_FORMAT "\n",
format);
return ok_format;
}
void Highs::clearZeroHessian() {
HighsHessian& hessian = model_.hessian_;
if (hessian.dim_) {
if (hessian.numNz() == 0) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Hessian has dimension %" HIGHSINT_FORMAT
" but no nonzeros, so is ignored\n",
hessian.dim_);
hessian.clear();
}
}
}
HighsStatus Highs::checkOptimality(const std::string& solver_type) {
assert(model_status_ == HighsModelStatus::kOptimal);
if (info_.num_primal_infeasibilities == 0 &&
info_.num_dual_infeasibilities <= 0)
return HighsStatus::kOk;
model_status_ = HighsModelStatus::kSolveError;
std::stringstream ss;
ss.str(std::string());
ss << highsFormatToString(
"%s solver claims optimality, but with num/max/sum "
"primal(%d/%g/%g)",
solver_type.c_str(), int(info_.num_primal_infeasibilities),
info_.max_primal_infeasibility, info_.sum_primal_infeasibilities);
if (info_.num_dual_infeasibilities > 0)
ss << highsFormatToString(
"and dual(%d/%g/%g)", int(info_.num_dual_infeasibilities),
info_.max_dual_infeasibility, info_.sum_dual_infeasibilities);
ss << " infeasibilities\n";
const std::string report_string = ss.str();
highsLogUser(options_.log_options, HighsLogType::kError, "%s",
report_string.c_str());
highsLogUser(options_.log_options, HighsLogType::kError,
"Setting model status to %s\n",
modelStatusToString(model_status_).c_str());
return HighsStatus::kError;
}
void Highs::callLpKktCheck(const HighsLp& lp, const std::string& message) {
lpKktCheck(this->model_status_, this->info_, lp, this->solution_,
this->basis_, this->options_, message);
}
HighsStatus Highs::invertRequirementError(std::string method_name) const {
assert(!ekk_instance_.status_.has_invert);
highsLogUser(options_.log_options, HighsLogType::kError,
"No invertible representation for %s\n", method_name.c_str());
return HighsStatus::kError;
}
HighsStatus Highs::handleInfCost() {
HighsLp& lp = this->model_.lp_;
if (!lp.has_infinite_cost_) return HighsStatus::kOk;
HighsLpMods& mods = lp.mods_;
double inf_cost = this->options_.infinite_cost;
for (HighsInt k = 0; k < 2; k++) {
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
double cost = lp.col_cost_[iCol];
if (cost > -inf_cost && cost < inf_cost) continue;
double lower = lp.col_lower_[iCol];
double upper = lp.col_upper_[iCol];
if (lp.isMip()) {
if (lp.integrality_[iCol] == HighsVarType::kInteger) {
lower = std::ceil(lower);
upper = std::floor(upper);
}
}
if (cost <= -inf_cost) {
if (lp.sense_ == ObjSense::kMinimize) {
if (upper < kHighsInf) {
if (k) lp.col_lower_[iCol] = upper;
} else {
highsLogUser(options_.log_options, HighsLogType::kError,
"Cannot minimize with a cost on variable %d of %g and "
"upper bound of %g\n",
int(iCol), cost, upper);
return HighsStatus::kError;
}
} else {
if (lower > -kHighsInf) {
if (k) lp.col_upper_[iCol] = lower;
} else {
highsLogUser(options_.log_options, HighsLogType::kError,
"Cannot maximize with a cost on variable %d of %g and "
"lower bound of %g\n",
int(iCol), cost, lower);
return HighsStatus::kError;
}
}
} else {
if (lp.sense_ == ObjSense::kMinimize) {
if (lower > -kHighsInf) {
if (k) lp.col_upper_[iCol] = lower;
} else {
highsLogUser(options_.log_options, HighsLogType::kError,
"Cannot minimize with a cost on variable %d of %g and "
"lower bound of %g\n",
int(iCol), cost, lower);
return HighsStatus::kError;
}
} else {
if (upper < kHighsInf) {
if (k) lp.col_lower_[iCol] = upper;
} else {
highsLogUser(options_.log_options, HighsLogType::kError,
"Cannot maximize with a cost on variable %d of %g and "
"upper bound of %g\n",
int(iCol), cost, upper);
return HighsStatus::kError;
}
}
}
if (k) {
mods.save_inf_cost_variable_index.push_back(iCol);
mods.save_inf_cost_variable_cost.push_back(cost);
mods.save_inf_cost_variable_lower.push_back(lower);
mods.save_inf_cost_variable_upper.push_back(upper);
lp.col_cost_[iCol] = 0;
}
}
}
lp.has_infinite_cost_ = false;
return HighsStatus::kOk;
}
HighsStatus Highs::optionChangeAction() {
if (this->iis_.valid_ && options_.iis_strategy != this->iis_.strategy_)
this->iis_.clear();
return HighsStatus::kOk;
}
void Highs::restoreInfCost(HighsStatus& return_status) {
HighsLp& lp = this->model_.lp_;
HighsBasis& basis = this->basis_;
HighsLpMods& mods = lp.mods_;
HighsInt num_inf_cost = mods.save_inf_cost_variable_index.size();
if (num_inf_cost <= 0) return;
assert(num_inf_cost);
for (HighsInt ix = 0; ix < num_inf_cost; ix++) {
HighsInt iCol = mods.save_inf_cost_variable_index[ix];
double cost = mods.save_inf_cost_variable_cost[ix];
double lower = mods.save_inf_cost_variable_lower[ix];
double upper = mods.save_inf_cost_variable_upper[ix];
double value = solution_.value_valid ? solution_.col_value[iCol] : 0;
if (basis.valid) {
assert(basis.col_status[iCol] != HighsBasisStatus::kBasic);
if (lp.col_lower_[iCol] == lower) {
basis.col_status[iCol] = HighsBasisStatus::kLower;
} else {
basis.col_status[iCol] = HighsBasisStatus::kUpper;
}
}
assert(lp.col_cost_[iCol] == 0);
if (value) this->info_.objective_function_value += value * cost;
lp.col_cost_[iCol] = cost;
lp.col_lower_[iCol] = lower;
lp.col_upper_[iCol] = upper;
}
lp.has_infinite_cost_ = true;
if (this->model_status_ == HighsModelStatus::kInfeasible) {
this->model_status_ = HighsModelStatus::kUnknown;
setHighsModelStatusAndClearSolutionAndBasis(this->model_status_);
return_status = highsStatusFromHighsModelStatus(model_status_);
}
}
HighsStatus Highs::userScale(HighsUserScaleData& data) {
if (!options_.user_objective_scale && !options_.user_bound_scale)
return HighsStatus::kOk;
initialiseUserScaleData(this->options_, data);
if (this->userScaleModel(data) == HighsStatus::kError)
return HighsStatus::kError;
this->userScaleSolution(data);
data.applied = true;
return HighsStatus::kOk;
}
HighsStatus Highs::userUnscale(HighsUserScaleData& data) {
if (!data.applied) return HighsStatus::kOk;
HighsStatus status = HighsStatus::kOk;
data.user_objective_scale *= -1;
data.user_bound_scale *= -1;
HighsStatus unscale_status = this->userScaleModel(data);
if (unscale_status == HighsStatus::kError) {
highsLogUser(
this->options_.log_options, HighsLogType::kError,
"Unexpected error removing user scaling from the incumbent model\n");
assert(unscale_status != HighsStatus::kError);
}
const bool update_kkt = true;
unscale_status = this->userScaleSolution(data, update_kkt);
highsLogUser(this->options_.log_options, HighsLogType::kInfo,
"After solving the user-scaled model, the unscaled solution "
"has objective value %.12g\n",
this->info_.objective_function_value);
if (model_status_ == HighsModelStatus::kOptimal &&
unscale_status != HighsStatus::kOk) {
highsLogUser(
this->options_.log_options, HighsLogType::kWarning,
"User scaled problem solved to optimality, but unscaled solution "
"does not satisfy feasibility and optimality tolerances\n");
status = HighsStatus::kWarning;
}
return status;
}
HighsStatus Highs::userScaleModel(HighsUserScaleData& data) {
userScaleLp(this->model_.lp_, data, false);
userScaleHessian(this->model_.hessian_, data, false);
HighsStatus return_status = userScaleStatus(this->options_.log_options, data);
if (return_status == HighsStatus::kError) return HighsStatus::kError;
userScaleLp(this->model_.lp_, data);
userScaleHessian(this->model_.hessian_, data);
return return_status;
}
HighsStatus Highs::userScaleSolution(HighsUserScaleData& data,
bool update_kkt) {
HighsStatus return_status = HighsStatus::kOk;
if (!data.user_objective_scale && !data.user_bound_scale)
return HighsStatus::kOk;
double objective_scale_value = std::pow(2, data.user_objective_scale);
double bound_scale_value = std::pow(2, data.user_bound_scale);
const HighsLp& lp = this->model_.lp_;
const bool has_integrality = lp.integrality_.size();
if (info_.primal_solution_status != kSolutionStatusNone) {
if (data.user_bound_scale) {
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
if (has_integrality &&
lp.integrality_[iCol] != HighsVarType::kContinuous)
continue;
this->solution_.col_value[iCol] *= bound_scale_value;
}
for (HighsInt iRow = 0; iRow < lp.num_row_; iRow++)
this->solution_.row_value[iRow] *= bound_scale_value;
}
}
if (info_.dual_solution_status != kSolutionStatusNone) {
if (data.user_objective_scale) {
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++)
this->solution_.col_dual[iCol] *= objective_scale_value;
for (HighsInt iRow = 0; iRow < lp.num_row_; iRow++)
this->solution_.row_dual[iRow] *= objective_scale_value;
}
}
if (!update_kkt) return return_status;
double objective_function_value =
info_.objective_function_value - model_.lp_.offset_;
objective_function_value *= (bound_scale_value * objective_scale_value);
objective_function_value += model_.lp_.offset_;
info_.objective_function_value = objective_function_value;
getKktFailures(options_, model_, solution_, basis_, info_);
return reportKktFailures(model_.lp_, options_, info_,
"After removing user scaling")
? HighsStatus::kWarning
: return_status;
}
void HighsIllConditioning::clear() { this->record.clear(); }
HighsStatus Highs::computeIllConditioning(
HighsIllConditioning& ill_conditioning, const bool constraint,
const HighsInt method, const double ill_conditioning_bound) {
const double kZeroMultiplier = 1e-6;
ill_conditioning.clear();
HighsLp& incumbent_lp = this->model_.lp_;
Highs conditioning;
const bool dev_conditioning = false;
conditioning.setOptionValue("output_flag", false); std::vector<HighsInt> basic_var;
HighsLp& ill_conditioning_lp = conditioning.model_.lp_;
if (method == 0) {
formIllConditioningLp0(ill_conditioning_lp, basic_var, constraint);
} else {
formIllConditioningLp1(ill_conditioning_lp, basic_var, constraint,
ill_conditioning_bound);
}
assert(assessLp(ill_conditioning_lp, this->options_) == HighsStatus::kOk);
HighsStatus return_status = conditioning.run();
HighsModelStatus model_status = conditioning.getModelStatus();
const std::string type = constraint ? "Constraint" : "Column";
const bool failed =
return_status != HighsStatus::kOk ||
(method == 0 && model_status != HighsModelStatus::kOptimal) ||
(method == 1 && (model_status != HighsModelStatus::kOptimal &&
model_status != HighsModelStatus::kInfeasible));
if (failed) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"\n%s view ill-conditioning analysis has failed\n",
type.c_str());
return HighsStatus::kError;
}
if (method == 1 && model_status == HighsModelStatus::kInfeasible) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"\n%s view ill-conditioning bound of %g is insufficient for "
"analysis: try %g\n",
type.c_str(), ill_conditioning_bound,
1e1 * ill_conditioning_bound);
return HighsStatus::kOk;
}
if (dev_conditioning) conditioning.writeSolution("", 1);
HighsSolution& solution = conditioning.solution_;
double multiplier_norm = 0;
for (HighsInt iRow = 0; iRow < incumbent_lp.num_row_; iRow++)
multiplier_norm += std::fabs(solution.col_value[iRow]);
assert(multiplier_norm > 0);
const double ill_conditioning_measure =
(method == 0 ? conditioning.getInfo().objective_function_value
: solution.row_value[conditioning.getNumRow() - 1]) /
multiplier_norm;
highsLogUser(
options_.log_options, HighsLogType::kInfo,
"\n%s view ill-conditioning analysis: 1-norm distance of basis matrix "
"from singularity is estimated to be %g\n",
type.c_str(), ill_conditioning_measure);
std::vector<std::pair<double, HighsInt>> abs_list;
for (HighsInt iRow = 0; iRow < incumbent_lp.num_row_; iRow++) {
double abs_multiplier =
std::fabs(solution.col_value[iRow]) / multiplier_norm;
if (abs_multiplier <= kZeroMultiplier) continue;
abs_list.push_back(std::make_pair(abs_multiplier, iRow));
}
std::sort(abs_list.begin(), abs_list.end());
std::stringstream ss;
const bool has_row_names =
HighsInt(incumbent_lp.row_names_.size()) == incumbent_lp.num_row_;
const bool has_col_names =
HighsInt(incumbent_lp.col_names_.size()) == incumbent_lp.num_col_;
const double coefficient_zero_tolerance = 1e-8;
auto printCoefficient = [&](const double multiplier, const bool first) {
if (std::fabs(multiplier) < coefficient_zero_tolerance) {
ss << "+ 0";
} else if (std::fabs(multiplier - 1) < coefficient_zero_tolerance) {
std::string str = first ? "" : "+ ";
ss << str;
} else if (std::fabs(multiplier + 1) < coefficient_zero_tolerance) {
std::string str = first ? "-" : "- ";
ss << str;
} else if (multiplier < 0) {
std::string str = first ? "-" : "- ";
ss << str << -multiplier << " ";
} else {
std::string str = first ? "" : "+ ";
ss << str << multiplier << " ";
}
};
for (HighsInt iX = int(abs_list.size()) - 1; iX >= 0; iX--) {
HighsInt iRow = abs_list[iX].second;
HighsIllConditioningRecord record;
record.index = iRow;
record.multiplier = solution.col_value[iRow] / multiplier_norm;
ill_conditioning.record.push_back(record);
}
HighsSparseMatrix& incumbent_matrix = incumbent_lp.a_matrix_;
if (constraint) {
HighsInt num_nz;
std::vector<HighsInt> index(incumbent_lp.num_col_);
std::vector<double> value(incumbent_lp.num_col_);
HighsInt* p_index = index.data();
double* p_value = value.data();
for (HighsInt iX = 0; iX < HighsInt(ill_conditioning.record.size()); iX++) {
ss.str(std::string());
bool newline = false;
HighsInt iRow = ill_conditioning.record[iX].index;
double multiplier = ill_conditioning.record[iX].multiplier;
num_nz = 0;
incumbent_matrix.getRow(iRow, num_nz, p_index, p_value);
std::string row_name = has_row_names ? incumbent_lp.row_names_[iRow]
: "R" + std::to_string(iRow);
ss << "(Mu=" << multiplier << ")" << row_name << ": ";
const double lower = incumbent_lp.row_lower_[iRow];
const double upper = incumbent_lp.row_upper_[iRow];
if (lower > -kHighsInf && lower != upper)
ss << incumbent_lp.row_lower_[iRow] << " <= ";
for (HighsInt iEl = 0; iEl < num_nz; iEl++) {
if (newline) {
ss << " ";
newline = false;
}
HighsInt iCol = index[iEl];
printCoefficient(value[iEl], iEl == 0);
std::string col_name = has_col_names ? incumbent_lp.col_names_[iCol]
: "C" + std::to_string(iCol);
ss << col_name << " ";
HighsInt length_ss = ss.str().length();
if (length_ss > 72 && iEl < num_nz - 1) {
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s\n",
ss.str().c_str());
ss.str(std::string());
newline = true;
}
}
if (upper < kHighsInf) {
if (lower == upper) {
ss << "= " << upper;
} else {
ss << "<= " << upper;
}
}
if (ss.str().length())
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s\n",
ss.str().c_str());
}
} else {
for (const auto& rec : ill_conditioning.record) {
ss.str(std::string());
bool newline = false;
double multiplier = rec.multiplier;
HighsInt iCol = basic_var[rec.index];
if (iCol < incumbent_lp.num_col_) {
std::string col_name = has_col_names ? incumbent_lp.col_names_[iCol]
: "C" + std::to_string(iCol);
ss << "(Mu=" << multiplier << ")" << col_name << ": ";
for (HighsInt iEl = incumbent_matrix.start_[iCol];
iEl < incumbent_matrix.start_[iCol + 1]; iEl++) {
if (newline) {
ss << " ";
newline = false;
} else {
if (iEl > incumbent_matrix.start_[iCol]) ss << " | ";
}
HighsInt iRow = incumbent_matrix.index_[iEl];
printCoefficient(incumbent_matrix.value_[iEl], true);
std::string row_name = has_row_names ? incumbent_lp.row_names_[iRow]
: "R" + std::to_string(iRow);
ss << row_name;
HighsInt length_ss = ss.str().length();
if (length_ss > 72 && iEl < incumbent_matrix.start_[iCol + 1] - 1) {
ss << " | ";
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s\n",
ss.str().c_str());
ss.str(std::string());
newline = true;
}
}
} else {
HighsInt iRow = iCol - incumbent_lp.num_col_;
std::string col_name = has_row_names
? "Slack_" + incumbent_lp.row_names_[iRow]
: "Slack_R" + std::to_string(iRow);
ss << "(Mu=" << multiplier << ")" << col_name << ": ";
}
if (ss.str().length())
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s\n",
ss.str().c_str());
}
}
return HighsStatus::kOk;
}
void Highs::formIllConditioningLp0(HighsLp& ill_conditioning_lp,
std::vector<HighsInt>& basic_var,
const bool constraint) {
HighsLp& incumbent_lp = this->model_.lp_;
ill_conditioning_lp.num_row_ = incumbent_lp.num_row_ + 1;
for (HighsInt iRow = 0; iRow < incumbent_lp.num_row_; iRow++) {
ill_conditioning_lp.row_lower_.push_back(0);
ill_conditioning_lp.row_upper_.push_back(0);
}
ill_conditioning_lp.row_lower_.push_back(1);
ill_conditioning_lp.row_upper_.push_back(1);
HighsSparseMatrix& incumbent_matrix = incumbent_lp.a_matrix_;
incumbent_matrix.ensureColwise();
HighsSparseMatrix& ill_conditioning_matrix = ill_conditioning_lp.a_matrix_;
ill_conditioning_matrix.num_row_ = ill_conditioning_lp.num_row_;
const HighsInt ill_conditioning_lp_e_row = ill_conditioning_lp.num_row_ - 1;
for (HighsInt iCol = 0; iCol < incumbent_lp.num_col_; iCol++) {
if (this->basis_.col_status[iCol] != HighsBasisStatus::kBasic) continue;
basic_var.push_back(iCol);
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(-kHighsInf);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
for (HighsInt iEl = incumbent_matrix.start_[iCol];
iEl < incumbent_matrix.start_[iCol + 1]; iEl++) {
ill_conditioning_matrix.index_.push_back(incumbent_matrix.index_[iEl]);
ill_conditioning_matrix.value_.push_back(incumbent_matrix.value_[iEl]);
}
if (!constraint) {
ill_conditioning_matrix.index_.push_back(ill_conditioning_lp_e_row);
ill_conditioning_matrix.value_.push_back(1.0);
}
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
}
for (HighsInt iRow = 0; iRow < incumbent_lp.num_row_; iRow++) {
if (this->basis_.row_status[iRow] != HighsBasisStatus::kBasic) continue;
basic_var.push_back(incumbent_lp.num_col_ + iRow);
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(-kHighsInf);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(iRow);
ill_conditioning_matrix.value_.push_back(-1.0);
if (!constraint) {
ill_conditioning_matrix.index_.push_back(ill_conditioning_lp_e_row);
ill_conditioning_matrix.value_.push_back(1.0);
}
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
}
if (constraint) {
for (HighsInt iRow = 0; iRow < incumbent_lp.num_row_; iRow++) {
ill_conditioning_matrix.index_.push_back(iRow);
ill_conditioning_matrix.value_.push_back(1.0);
}
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_matrix.num_row_ = incumbent_lp.num_row_;
ill_conditioning_matrix.num_col_ = incumbent_lp.num_row_ + 1;
ill_conditioning_matrix.ensureRowwise();
ill_conditioning_matrix.format_ = MatrixFormat::kColwise;
}
for (HighsInt iRow = 0; iRow < incumbent_lp.num_row_; iRow++) {
ill_conditioning_lp.col_cost_.push_back(1);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(iRow);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.col_cost_.push_back(1);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(iRow);
ill_conditioning_matrix.value_.push_back(-1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
}
ill_conditioning_lp.num_col_ = 3 * incumbent_lp.num_row_;
ill_conditioning_matrix.num_col_ = ill_conditioning_lp.num_col_;
ill_conditioning_matrix.num_row_ = ill_conditioning_lp.num_row_;
}
void Highs::formIllConditioningLp1(HighsLp& ill_conditioning_lp,
std::vector<HighsInt>& basic_var,
const bool constraint,
const double ill_conditioning_bound) {
HighsLp& incumbent_lp = this->model_.lp_;
const HighsInt incumbent_num_row = incumbent_lp.num_row_;
const HighsInt c4_offset = 0;
const HighsInt c1_offset = incumbent_num_row;
const HighsInt c7_offset = 2 * incumbent_num_row;
const HighsInt c6_offset = 3 * incumbent_num_row;
const HighsInt c5_offset = 3 * incumbent_num_row + 1;
for (HighsInt iRow = 0; iRow < c6_offset; iRow++) {
ill_conditioning_lp.row_lower_.push_back(0);
ill_conditioning_lp.row_upper_.push_back(0);
}
HighsSparseMatrix& incumbent_matrix = incumbent_lp.a_matrix_;
incumbent_matrix.ensureColwise();
HighsSparseMatrix& ill_conditioning_matrix = ill_conditioning_lp.a_matrix_;
ill_conditioning_lp.num_col_ = 0;
for (HighsInt iCol = 0; iCol < incumbent_lp.num_col_; iCol++) {
if (this->basis_.col_status[iCol] != HighsBasisStatus::kBasic) continue;
basic_var.push_back(iCol);
ill_conditioning_lp.col_names_.push_back(
"y_" + std::to_string(ill_conditioning_lp.num_col_));
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(-kHighsInf);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
for (HighsInt iEl = incumbent_matrix.start_[iCol];
iEl < incumbent_matrix.start_[iCol + 1]; iEl++) {
ill_conditioning_matrix.index_.push_back(incumbent_matrix.index_[iEl]);
ill_conditioning_matrix.value_.push_back(incumbent_matrix.value_[iEl]);
}
if (!constraint) {
ill_conditioning_matrix.index_.push_back(c1_offset +
ill_conditioning_lp.num_col_);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.index_.push_back(c6_offset);
ill_conditioning_matrix.value_.push_back(1.0);
}
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
}
for (HighsInt iRow = 0; iRow < incumbent_num_row; iRow++) {
if (this->basis_.row_status[iRow] != HighsBasisStatus::kBasic) continue;
basic_var.push_back(incumbent_lp.num_col_ + iRow);
ill_conditioning_lp.col_names_.push_back(
"y_" + std::to_string(ill_conditioning_lp.num_col_));
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(-kHighsInf);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(iRow);
ill_conditioning_matrix.value_.push_back(-1.0);
if (!constraint) {
ill_conditioning_matrix.index_.push_back(c1_offset +
ill_conditioning_lp.num_col_);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.index_.push_back(c6_offset);
ill_conditioning_matrix.value_.push_back(1.0);
}
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
}
assert(ill_conditioning_lp.num_col_ == incumbent_num_row);
if (constraint) {
for (HighsInt iRow = 0; iRow < incumbent_num_row; iRow++) {
ill_conditioning_matrix.index_.push_back(iRow);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
}
for (HighsInt iRow = 0; iRow < incumbent_num_row; iRow++)
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
for (HighsInt iRow = 0; iRow < incumbent_num_row; iRow++) {
ill_conditioning_matrix.index_.push_back(iRow);
ill_conditioning_matrix.value_.push_back(1.0);
}
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_matrix.num_col_ = c6_offset + 1;
ill_conditioning_matrix.num_row_ = incumbent_num_row;
ill_conditioning_matrix.ensureRowwise();
ill_conditioning_matrix.format_ = MatrixFormat::kColwise;
ill_conditioning_matrix.num_col_ = incumbent_num_row;
ill_conditioning_matrix.num_row_ = c6_offset + 1;
}
assert(ill_conditioning_lp.num_col_ == incumbent_num_row);
ill_conditioning_lp.num_row_ = 3 * incumbent_num_row + 2;
for (HighsInt iRow = 0; iRow < incumbent_num_row; iRow++) {
ill_conditioning_lp.col_names_.push_back("u_" + std::to_string(iRow));
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(c1_offset + iRow);
ill_conditioning_matrix.value_.push_back(-1.0);
ill_conditioning_matrix.index_.push_back(c7_offset + iRow);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
ill_conditioning_lp.col_names_.push_back("w_" + std::to_string(iRow));
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(c1_offset + iRow);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.index_.push_back(c7_offset + iRow);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
}
for (HighsInt iRow = 0; iRow < incumbent_num_row; iRow++) {
ill_conditioning_lp.col_names_.push_back("s_" + std::to_string(iRow));
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(c4_offset + iRow);
ill_conditioning_matrix.value_.push_back(-1.0);
ill_conditioning_matrix.index_.push_back(c5_offset);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
ill_conditioning_lp.col_names_.push_back("t_" + std::to_string(iRow));
ill_conditioning_lp.col_cost_.push_back(0);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(c4_offset + iRow);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.index_.push_back(c5_offset);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
}
ill_conditioning_lp.row_lower_.push_back(1);
ill_conditioning_lp.row_upper_.push_back(1);
assert(ill_conditioning_bound > 0);
ill_conditioning_lp.row_lower_.push_back(-kHighsInf);
ill_conditioning_lp.row_upper_.push_back(ill_conditioning_bound);
assert(HighsInt(ill_conditioning_lp.row_lower_.size()) ==
ill_conditioning_lp.num_row_);
assert(HighsInt(ill_conditioning_lp.row_upper_.size()) ==
ill_conditioning_lp.num_row_);
for (HighsInt iRow = 0; iRow < incumbent_num_row; iRow++) {
ill_conditioning_lp.col_names_.push_back("IfsPlus_" + std::to_string(iRow));
ill_conditioning_lp.col_cost_.push_back(1);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(c7_offset + iRow);
ill_conditioning_matrix.value_.push_back(-1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
ill_conditioning_lp.col_names_.push_back("IfsMinus_" +
std::to_string(iRow));
ill_conditioning_lp.col_cost_.push_back(1);
ill_conditioning_lp.col_lower_.push_back(0);
ill_conditioning_lp.col_upper_.push_back(kHighsInf);
ill_conditioning_matrix.index_.push_back(c7_offset + iRow);
ill_conditioning_matrix.value_.push_back(1.0);
ill_conditioning_matrix.start_.push_back(
HighsInt(ill_conditioning_matrix.index_.size()));
ill_conditioning_lp.num_col_++;
}
assert(ill_conditioning_lp.num_col_ == 7 * incumbent_num_row);
assert(ill_conditioning_lp.num_row_ == 3 * incumbent_num_row + 2);
ill_conditioning_matrix.num_col_ = ill_conditioning_lp.num_col_;
ill_conditioning_matrix.num_row_ = ill_conditioning_lp.num_row_;
}
bool Highs::infeasibleBoundsOk() {
const HighsLogOptions& log_options = this->options_.log_options;
HighsLp& lp = this->model_.lp_;
HighsInt num_true_infeasible_bound = 0;
HighsInt num_ok_infeasible_bound = 0;
const bool has_integrality = lp.integrality_.size() > 0;
bool performed_inward_integer_rounding = false;
auto infeasibleBoundOk = [&](const std::string& type, const HighsInt iX,
double& lower, double& upper) {
double range = upper - lower;
assert(range < 0);
if (range >= 0) return true;
if (range > -this->options_.primal_feasibility_tolerance) {
num_ok_infeasible_bound++;
bool report = num_ok_infeasible_bound <= 10;
bool integer_lower = lower == std::floor(lower + 0.5);
bool integer_upper = upper == std::floor(upper + 0.5);
assert(!integer_lower || !integer_upper);
if (integer_lower) {
if (report)
highsLogUser(log_options, HighsLogType::kInfo,
"%s %d bounds [%g, %g] have infeasibility = %g so set "
"upper bound to %g\n",
type.c_str(), int(iX), lower, upper, range, lower);
upper = lower;
} else if (integer_upper) {
if (report)
highsLogUser(log_options, HighsLogType::kInfo,
"%s %d bounds [%g, %g] have infeasibility = %g so set "
"lower bound to %g\n",
type.c_str(), int(iX), lower, upper, range, upper);
lower = upper;
} else {
double mid = 0.5 * (lower + upper);
if (report)
highsLogUser(log_options, HighsLogType::kInfo,
"%s %d bounds [%g, %g] have infeasibility = %g so set "
"both bounds to %g\n",
type.c_str(), int(iX), lower, upper, range, mid);
lower = mid;
upper = mid;
}
return true;
}
num_true_infeasible_bound++;
if (num_true_infeasible_bound <= 10)
highsLogUser(
log_options, HighsLogType::kInfo,
"%s %d bounds [%g, %g] have excessive infeasibility = %g%s\n",
type.c_str(), int(iX), lower, upper, range,
performed_inward_integer_rounding ? " due to inward integer rounding"
: "");
return false;
};
const bool perform_inward_integer_rounding = !this->options_.solve_relaxation;
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
performed_inward_integer_rounding = false;
double lower = lp.col_lower_[iCol];
double upper = lp.col_upper_[iCol];
if (has_integrality) {
if (lp.integrality_[iCol] == HighsVarType::kSemiContinuous ||
lp.integrality_[iCol] == HighsVarType::kSemiInteger)
continue;
if (perform_inward_integer_rounding &&
lp.integrality_[iCol] == HighsVarType::kInteger) {
double feastol = this->options_.mip_feasibility_tolerance;
double integer_lower = std::ceil(lower - feastol);
double integer_upper = std::floor(upper + feastol);
assert(integer_lower >= lower);
assert(integer_upper <= upper);
performed_inward_integer_rounding =
integer_lower > lower || integer_upper < upper;
lower = integer_lower;
upper = integer_upper;
}
}
if (lower > upper) {
if (infeasibleBoundOk("Column", iCol, lower, upper)) {
lp.col_lower_[iCol] = lower;
lp.col_upper_[iCol] = upper;
}
}
}
performed_inward_integer_rounding = false;
for (HighsInt iRow = 0; iRow < lp.num_row_; iRow++) {
if (lp.row_lower_[iRow] > lp.row_upper_[iRow])
infeasibleBoundOk("Row", iRow, lp.row_lower_[iRow], lp.row_upper_[iRow]);
}
if (num_ok_infeasible_bound > 0)
highsLogUser(log_options, HighsLogType::kInfo,
"Model has %d small inconsistent bound(s): rectified\n",
int(num_ok_infeasible_bound));
if (num_true_infeasible_bound > 0)
highsLogUser(log_options, HighsLogType::kInfo,
"Model has %d significant inconsistent bound(s): infeasible\n",
int(num_true_infeasible_bound));
return num_true_infeasible_bound == 0;
}
bool Highs::validLinearObjective(const HighsLinearObjective& linear_objective,
const HighsInt iObj) const {
HighsInt linear_objective_coefficients_size =
linear_objective.coefficients.size();
if (linear_objective_coefficients_size != this->model_.lp_.num_col_) {
highsLogUser(
options_.log_options, HighsLogType::kError,
"Coefficient vector for linear objective %s has size %d != %d = "
"lp.num_col_\n",
iObj >= 0 ? std::to_string(iObj).c_str() : "",
int(linear_objective_coefficients_size),
int(this->model_.lp_.num_col_));
return false;
}
if (!options_.blend_multi_objectives &&
hasRepeatedLinearObjectivePriorities(&linear_objective)) {
highsLogUser(
options_.log_options, HighsLogType::kError,
"Repeated priorities for lexicographic optimization is illegal\n");
return false;
}
return true;
}
bool Highs::hasRepeatedLinearObjectivePriorities(
const HighsLinearObjective* linear_objective) const {
HighsInt num_linear_objective = this->multi_linear_objective_.size();
if (num_linear_objective <= 0 ||
(num_linear_objective <= 1 && !linear_objective))
return false;
for (HighsInt iObj0 = 0; iObj0 < num_linear_objective; iObj0++) {
HighsInt priority0 = this->multi_linear_objective_[iObj0].priority;
for (HighsInt iObj1 = iObj0 + 1; iObj1 < num_linear_objective; iObj1++) {
HighsInt priority1 = this->multi_linear_objective_[iObj1].priority;
if (priority1 == priority0) return true;
}
if (linear_objective) {
if (linear_objective->priority == priority0) return true;
}
}
return false;
}
static bool comparison(std::pair<HighsInt, HighsInt> x1,
std::pair<HighsInt, HighsInt> x2) {
return x1.first > x2.first;
}
HighsStatus Highs::returnFromLexicographicOptimization(
HighsStatus return_status, HighsInt original_lp_num_row) {
HighsModelStatus model_status = this->model_status_;
HighsInfo info = this->info_;
HighsInt num_linear_objective = this->multi_linear_objective_.size();
if (num_linear_objective > 1) {
this->deleteRows(original_lp_num_row, this->model_.lp_.num_row_ - 1);
this->model_status_ = model_status;
this->info_ = info;
info_.objective_function_value = 0;
info_.basis_validity = kBasisValidityInvalid;
info_.invalidateDualKkt();
this->solution_.value_valid = true;
this->model_.lp_.col_cost_.assign(this->model_.lp_.num_col_, 0);
}
return return_status;
}
HighsStatus Highs::multiobjectiveSolve() {
const HighsInt coeff_logging_size_limit = 10;
HighsInt num_linear_objective = this->multi_linear_objective_.size();
assert(num_linear_objective > 0);
HighsLp& lp = this->model_.lp_;
for (HighsInt iObj = 0; iObj < num_linear_objective; iObj++) {
HighsLinearObjective& multi_linear_objective =
this->multi_linear_objective_[iObj];
if (multi_linear_objective.coefficients.size() !=
static_cast<size_t>(lp.num_col_)) {
highsLogUser(options_.log_options, HighsLogType::kError,
"Multiple linear objective coefficient vector %d has size "
"incompatible with model\n",
int(iObj));
return HighsStatus::kError;
}
}
std::unique_ptr<std::stringstream> multi_objective_log;
highsLogUser(options_.log_options, HighsLogType::kInfo,
"Solving with %d multiple linear objectives, %s\n",
int(num_linear_objective),
this->options_.blend_multi_objectives
? "blending objectives by weight"
: "using lexicographic optimization by priority");
highsLogUser(
options_.log_options, HighsLogType::kInfo,
"Ix weight offset abs_tol rel_tol priority%s\n",
lp.num_col_ < coeff_logging_size_limit ? " coefficients" : "");
for (HighsInt iObj = 0; iObj < num_linear_objective; iObj++) {
HighsLinearObjective& linear_objective =
this->multi_linear_objective_[iObj];
multi_objective_log =
std::unique_ptr<std::stringstream>(new std::stringstream());
*multi_objective_log << highsFormatToString(
"%2d %11.6g %11.6g %11.6g %11.6g %11d ", int(iObj),
linear_objective.weight, linear_objective.offset,
linear_objective.abs_tolerance, linear_objective.rel_tolerance,
int(linear_objective.priority));
if (lp.num_col_ < coeff_logging_size_limit) {
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++)
*multi_objective_log << highsFormatToString(
"%s c_{%1d} = %g", iCol == 0 ? "" : ",", int(iCol),
linear_objective.coefficients[iCol]);
}
*multi_objective_log << "\n";
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s",
multi_objective_log->str().c_str());
}
this->clearSolverDualData();
if (this->options_.blend_multi_objectives) {
lp.offset_ = 0;
lp.col_cost_.assign(lp.num_col_, 0);
for (HighsInt iObj = 0; iObj < num_linear_objective; iObj++) {
HighsLinearObjective& multi_linear_objective =
this->multi_linear_objective_[iObj];
lp.offset_ +=
multi_linear_objective.weight * multi_linear_objective.offset;
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++)
lp.col_cost_[iCol] += multi_linear_objective.weight *
multi_linear_objective.coefficients[iCol];
}
lp.sense_ = ObjSense::kMinimize;
multi_objective_log =
std::unique_ptr<std::stringstream>(new std::stringstream());
*multi_objective_log << highsFormatToString(
"Solving with blended objective");
if (lp.num_col_ < coeff_logging_size_limit) {
*multi_objective_log << highsFormatToString(
": %s %g", lp.sense_ == ObjSense::kMinimize ? "min" : "max",
lp.offset_);
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++)
*multi_objective_log << highsFormatToString(
" + (%g) x[%d]", lp.col_cost_[iCol], int(iCol));
}
*multi_objective_log << "\n";
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s",
multi_objective_log->str().c_str());
return this->optimizeModel();
}
if (model_.isQp() && num_linear_objective > 1) {
highsLogUser(
options_.log_options, HighsLogType::kError,
"Cannot perform non-trivial lexicographic optimization for QP\n");
return HighsStatus::kError;
}
if (hasRepeatedLinearObjectivePriorities()) {
highsLogUser(
options_.log_options, HighsLogType::kError,
"Repeated priorities for lexicographic optimization is illegal\n");
return HighsStatus::kError;
}
std::vector<std::pair<HighsInt, HighsInt>> priority_objective;
for (HighsInt iObj = 0; iObj < num_linear_objective; iObj++)
priority_objective.push_back(
std::make_pair(this->multi_linear_objective_[iObj].priority, iObj));
std::sort(priority_objective.begin(), priority_objective.end(), comparison);
lp.offset_ = 0;
lp.col_cost_.assign(lp.num_col_, 0);
const HighsInt original_lp_num_row = lp.num_row_;
std::vector<HighsInt> index(lp.num_col_);
std::vector<double> value(lp.num_col_);
HighsSolution solution;
for (HighsInt iIx = 0; iIx < num_linear_objective; iIx++) {
HighsInt priority = priority_objective[iIx].first;
HighsInt iObj = priority_objective[iIx].second;
HighsLinearObjective& linear_objective =
this->multi_linear_objective_[iObj];
lp.offset_ = linear_objective.offset;
lp.col_cost_ = linear_objective.coefficients;
lp.sense_ =
linear_objective.weight > 0 ? ObjSense::kMinimize : ObjSense::kMaximize;
if (lp.isMip() && solution.value_valid) {
HighsStatus set_solution_status = this->setSolution(solution);
if (set_solution_status == HighsStatus::kError) {
highsLogUser(options_.log_options, HighsLogType::kError,
"Failure to use one MIP to provide an integer feasible "
"solution of the next\n");
return returnFromLexicographicOptimization(HighsStatus::kError,
original_lp_num_row);
}
bool valid, integral, feasible;
HighsStatus assess_primal_solution =
assessPrimalSolution(valid, integral, feasible);
if (!valid || !integral || !feasible) {
highsLogUser(options_.log_options, HighsLogType::kWarning,
"Failure to use one MIP to provide an integer feasible "
"solution of the next: "
"status is valid = %s, integral = %s, feasible = %s\n",
highsBoolToString(valid).c_str(),
highsBoolToString(integral).c_str(),
highsBoolToString(feasible).c_str());
}
}
multi_objective_log =
std::unique_ptr<std::stringstream>(new std::stringstream());
*multi_objective_log << highsFormatToString("Solving with objective %d",
int(iObj));
if (lp.num_col_ < coeff_logging_size_limit) {
*multi_objective_log << highsFormatToString(
": %s %g", lp.sense_ == ObjSense::kMinimize ? "min" : "max",
lp.offset_);
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++)
*multi_objective_log << highsFormatToString(
" + (%g) x[%d]", lp.col_cost_[iCol], int(iCol));
}
*multi_objective_log << "\n";
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s",
multi_objective_log->str().c_str());
HighsStatus optimize_model_status = this->optimizeModel();
if (optimize_model_status == HighsStatus::kError)
return returnFromLexicographicOptimization(HighsStatus::kError,
original_lp_num_row);
if (model_status_ != HighsModelStatus::kOptimal) {
highsLogUser(options_.log_options, HighsLogType::kWarning,
"After priority %d solve, model status is %s\n",
int(priority), modelStatusToString(model_status_).c_str());
return returnFromLexicographicOptimization(HighsStatus::kWarning,
original_lp_num_row);
}
if (iIx == num_linear_objective - 1) break;
if (lp.isMip()) {
solution.col_value = this->solution_.col_value;
solution.value_valid = true;
}
HighsInt nnz = 0;
for (HighsInt iCol = 0; iCol < lp.num_col_; iCol++) {
if (lp.col_cost_[iCol]) {
index[nnz] = iCol;
value[nnz] = lp.col_cost_[iCol];
nnz++;
}
}
double objective = info_.objective_function_value;
HighsStatus add_row_status;
double lower_bound = -kHighsInf;
double upper_bound = kHighsInf;
if (lp.sense_ == ObjSense::kMinimize) {
if (linear_objective.abs_tolerance >= 0)
upper_bound = objective + linear_objective.abs_tolerance;
if (linear_objective.rel_tolerance >= 0) {
if (objective >= 0) {
upper_bound = std::min(
objective * (1.0 + linear_objective.rel_tolerance), upper_bound);
} else if (objective < 0) {
upper_bound = std::min(
objective * (1.0 - linear_objective.rel_tolerance), upper_bound);
}
}
upper_bound -= lp.offset_;
} else {
if (linear_objective.abs_tolerance >= 0)
lower_bound = objective - linear_objective.abs_tolerance;
if (linear_objective.rel_tolerance >= 0) {
if (objective >= 0) {
lower_bound = std::max(
objective * (1.0 - linear_objective.rel_tolerance), lower_bound);
} else if (objective < 0) {
lower_bound = std::max(
objective * (1.0 + linear_objective.rel_tolerance), lower_bound);
}
}
lower_bound -= lp.offset_;
}
if (lower_bound == -kHighsInf && upper_bound == kHighsInf)
highsLogUser(options_.log_options, HighsLogType::kWarning,
"After priority %d solve, no objective constraint due to "
"absolute tolerance being %g < 0,"
" and relative tolerance being %g < 0\n",
int(priority), linear_objective.abs_tolerance,
linear_objective.rel_tolerance);
multi_objective_log =
std::unique_ptr<std::stringstream>(new std::stringstream());
*multi_objective_log << highsFormatToString(
"Add constraint for objective %d: ", int(iObj));
if (nnz < coeff_logging_size_limit) {
*multi_objective_log << highsFormatToString("%g <= ", lower_bound);
for (HighsInt iEl = 0; iEl < nnz; iEl++)
*multi_objective_log << highsFormatToString(
"%s(%g) x[%d]", iEl > 0 ? " + " : "", value[iEl], int(index[iEl]));
*multi_objective_log << highsFormatToString(" <= %g\n", upper_bound);
} else {
*multi_objective_log << highsFormatToString("Bounds [%g, %g]\n",
lower_bound, upper_bound);
}
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s",
multi_objective_log->str().c_str());
add_row_status =
this->addRow(lower_bound, upper_bound, nnz, index.data(), value.data());
assert(add_row_status == HighsStatus::kOk);
}
return returnFromLexicographicOptimization(HighsStatus::kOk,
original_lp_num_row);
}
bool Highs::tryPdlpCleanup(HighsInt& pdlp_cleanup_iteration_limit,
const HighsInfo& presolved_lp_info) const {
const double tolerance_margin = 1e2;
bool no_cleanup = false;
double max_relative_violation = 0;
auto noCleanup = [&](const std::string& kkt_name, const double kkt_error,
const double kkt_tolerance) {
double use_kkt_tolerance =
this->options_.kkt_tolerance != kDefaultKktTolerance
? this->options_.kkt_tolerance
: kkt_tolerance;
double relative_violation = kkt_error / use_kkt_tolerance;
if (relative_violation > tolerance_margin)
printf(
"KKT measure (%11.4g, %11.4g) gives relative violation of %11.4g for "
"%s\n",
kkt_error, use_kkt_tolerance, relative_violation, kkt_name.c_str());
max_relative_violation =
std::max(relative_violation, max_relative_violation);
no_cleanup = max_relative_violation > tolerance_margin;
};
noCleanup("Max relative primal infeasibility",
this->info_.max_relative_primal_infeasibility,
this->options_.primal_feasibility_tolerance);
noCleanup("Max relative dual infeasibility",
this->info_.max_relative_dual_infeasibility,
this->options_.dual_feasibility_tolerance);
noCleanup("Max relative primal residual error",
this->info_.max_relative_primal_residual_error,
this->options_.primal_residual_tolerance);
noCleanup("Max relative dual residual error",
this->info_.max_relative_dual_residual_error,
this->options_.dual_residual_tolerance);
noCleanup("Primal-dual objective error",
this->info_.primal_dual_objective_error,
this->options_.optimality_tolerance);
if (no_cleanup) {
highsLogUser(options_.log_options, HighsLogType::kInfo,
"No PDLP cleanup due to KKT errors exceeding tolerances by a "
"max factor = %g > %g = allowed margin\n",
max_relative_violation, tolerance_margin);
return false;
}
if (presolved_lp_info.pdlp_iteration_count > 0) {
HighsInt ten_percent_pdlp_iteration_count =
presolved_lp_info.pdlp_iteration_count / 10;
pdlp_cleanup_iteration_limit =
std::max(HighsInt(10000), ten_percent_pdlp_iteration_count);
} else {
pdlp_cleanup_iteration_limit = 1000;
}
return true;
}
void HighsLinearObjective::clear() {
this->weight = 0.0;
this->offset = 0.0;
this->coefficients.clear();
this->abs_tolerance = 0.0;
this->rel_tolerance = 0.0;
this->priority = 0;
}
void HighsSubSolverCallTime::initialise() {
this->num_call.assign(kSubSolverCount, 0);
this->run_time.assign(kSubSolverCount, 0);
this->name.assign(kSubSolverCount, "");
this->name[kSubSolverDuSimplexBasis] = "Du simplex (basis)";
this->name[kSubSolverDuSimplexNoBasis] = "Du simplex (no basis)";
this->name[kSubSolverPrSimplexBasis] = "Pr simplex (basis)";
this->name[kSubSolverPrSimplexNoBasis] = "Pr simplex (no basis)";
this->name[kSubSolverHipo] = "HiPO";
this->name[kSubSolverIpx] = "IPX";
this->name[kSubSolverHipoAc] = "HiPO (AC)";
this->name[kSubSolverIpxAc] = "IPX (AC)";
this->name[kSubSolverPdlp] = "PDLP";
this->name[kSubSolverQpAsm] = "QP ASM";
this->name[kSubSolverMip] = "MIP";
this->name[kSubSolverSubMip] = "Sub-MIP";
}
void HighsSubSolverCallTime::add(
const HighsSubSolverCallTime& sub_solver_call_time,
const bool analytic_centre) {
for (HighsInt Ix = 0; Ix < kSubSolverCount; Ix++) {
HighsInt ToIx = Ix;
if (Ix == kSubSolverHipo) {
if (analytic_centre) ToIx = kSubSolverHipoAc;
} else if (Ix == kSubSolverIpx) {
if (analytic_centre) ToIx = kSubSolverIpxAc;
}
this->num_call[ToIx] += sub_solver_call_time.num_call[Ix];
this->run_time[ToIx] += sub_solver_call_time.run_time[Ix];
}
}
void Highs::reportSubSolverCallTime() const {
double mip_time = this->sub_solver_call_time_.run_time[kSubSolverMip];
std::stringstream ss;
ss.str(std::string());
ss << highsFormatToString(
"\nSub-solver timing\nSolver Calls Time "
"Time/call");
if (mip_time > 0) ss << " MIP%";
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s\n",
ss.str().c_str());
double sum_mip_sub_solve_time = 0;
for (HighsInt Ix = 0; Ix < kSubSolverCount; Ix++) {
if (this->sub_solver_call_time_.num_call[Ix]) {
ss.str(std::string());
ss << highsFormatToString(
"%-21s %9d %11.4e %11.4e",
this->sub_solver_call_time_.name[Ix].c_str(),
int(this->sub_solver_call_time_.num_call[Ix]),
this->sub_solver_call_time_.run_time[Ix],
this->sub_solver_call_time_.run_time[Ix] /
(1.0 * this->sub_solver_call_time_.num_call[Ix]));
if (mip_time > 0 && Ix != kSubSolverMip) {
if (Ix != kSubSolverHipoAc && Ix != kSubSolverIpxAc)
sum_mip_sub_solve_time += this->sub_solver_call_time_.run_time[Ix];
ss << highsFormatToString(
" %5.1f",
1e2 * this->sub_solver_call_time_.run_time[Ix] / mip_time);
}
highsLogUser(options_.log_options, HighsLogType::kInfo, "%s\n",
ss.str().c_str());
}
}
if (mip_time > 0)
highsLogUser(options_.log_options, HighsLogType::kInfo,
"TOTAL (excluding AC) %11.4e %5.1f\n",
sum_mip_sub_solve_time,
1e2 * sum_mip_sub_solve_time / mip_time);
}