Google OR-Tools v9.11
a fast and portable software suite for combinatorial optimization
Loading...
Searching...
No Matches
quadratic_program.cc
Go to the documentation of this file.
1// Copyright 2010-2024 Google LLC
2// Licensed under the Apache License, Version 2.0 (the "License");
3// you may not use this file except in compliance with the License.
4// You may obtain a copy of the License at
5//
6// http://www.apache.org/licenses/LICENSE-2.0
7//
8// Unless required by applicable law or agreed to in writing, software
9// distributed under the License is distributed on an "AS IS" BASIS,
10// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
11// See the License for the specific language governing permissions and
12// limitations under the License.
13
15
16#include <algorithm>
17#include <cstdint>
18#include <limits>
19#include <optional>
20#include <string>
21#include <tuple>
22#include <vector>
23
24#include "Eigen/Core"
25#include "Eigen/SparseCore"
26#include "absl/log/check.h"
27#include "absl/status/status.h"
28#include "absl/status/statusor.h"
29#include "absl/strings/str_cat.h"
30#include "absl/strings/string_view.h"
32#include "ortools/linear_solver/linear_solver.pb.h"
33
35
36using ::Eigen::VectorXd;
37
39 const int64_t var_lb_size = qp.variable_lower_bounds.size();
40 const int64_t con_lb_size = qp.constraint_lower_bounds.size();
41
42 if (var_lb_size != qp.variable_upper_bounds.size()) {
43 return absl::InvalidArgumentError(absl::StrCat(
44 "Inconsistent dimensions: variable lower bound vector has size ",
45 var_lb_size, " while variable upper bound vector has size ",
46 qp.variable_upper_bounds.size()));
47 }
48 if (var_lb_size != qp.objective_vector.size()) {
49 return absl::InvalidArgumentError(absl::StrCat(
50 "Inconsistent dimensions: variable lower bound vector has size ",
51 var_lb_size, " while objective vector has size ",
52 qp.objective_vector.size()));
53 }
54 if (var_lb_size != qp.constraint_matrix.cols()) {
55 return absl::InvalidArgumentError(absl::StrCat(
56 "Inconsistent dimensions: variable lower bound vector has size ",
57 var_lb_size, " while constraint matrix has ",
58 qp.constraint_matrix.cols(), " columns"));
59 }
60 if (qp.objective_matrix.has_value() &&
61 var_lb_size != qp.objective_matrix->rows()) {
62 return absl::InvalidArgumentError(absl::StrCat(
63 "Inconsistent dimensions: variable lower bound vector has size ",
64 var_lb_size, " while objective matrix has ",
65 qp.objective_matrix->rows(), " rows"));
66 }
67 if (con_lb_size != qp.constraint_upper_bounds.size()) {
68 return absl::InvalidArgumentError(absl::StrCat(
69 "Inconsistent dimensions: constraint lower bound vector has size ",
70 con_lb_size, " while constraint upper bound vector has size ",
71 qp.constraint_upper_bounds.size()));
72 }
73 if (con_lb_size != qp.constraint_matrix.rows()) {
74 return absl::InvalidArgumentError(absl::StrCat(
75 "Inconsistent dimensions: constraint lower bound vector has size ",
76 con_lb_size, " while constraint matrix has ",
77 qp.constraint_matrix.rows(), " rows "));
78 }
79
80 if (qp.variable_names.has_value() &&
81 var_lb_size != qp.variable_names->size()) {
82 return absl::InvalidArgumentError(absl::StrCat(
83 "Inconsistent dimensions: variable lower bound vector has size ",
84 var_lb_size, " while variable names has size ",
85 qp.variable_names->size()));
86 }
87
88 if (qp.constraint_names.has_value() &&
89 con_lb_size != qp.constraint_names->size()) {
90 return absl::InvalidArgumentError(absl::StrCat(
91 "Inconsistent dimensions: constraint lower bound vector has size ",
92 con_lb_size, " while constraint names has size ",
93 qp.constraint_names->size()));
94 }
95
96 return absl::OkStatus();
97}
98
99absl::StatusOr<QuadraticProgram> QpFromMpModelProto(
100 const MPModelProto& proto, bool relax_integer_variables,
101 bool include_names) {
102 if (!proto.general_constraint().empty()) {
103 return absl::InvalidArgumentError("General constraints are not supported.");
104 }
105 const int primal_size = proto.variable_size();
106 const int dual_size = proto.constraint_size();
107 QuadraticProgram qp(primal_size, dual_size);
108 if (include_names) {
109 qp.problem_name = proto.name();
110 qp.variable_names = std::vector<std::string>(primal_size);
111 qp.constraint_names = std::vector<std::string>(dual_size);
112 }
113 for (int i = 0; i < primal_size; ++i) {
114 const auto& var = proto.variable(i);
115 qp.variable_lower_bounds[i] = var.lower_bound();
116 qp.variable_upper_bounds[i] = var.upper_bound();
117 qp.objective_vector[i] = var.objective_coefficient();
118 if (var.is_integer() && !relax_integer_variables) {
119 return absl::InvalidArgumentError(
120 "Integer variable encountered with relax_integer_variables == false");
121 }
122 if (include_names) {
123 (*qp.variable_names)[i] = var.name();
124 }
125 }
126 std::vector<int> nonzeros_by_column(primal_size);
127 for (int i = 0; i < dual_size; ++i) {
128 const auto& con = proto.constraint(i);
129 for (int j = 0; j < con.var_index_size(); ++j) {
130 if (con.var_index(j) < 0 || con.var_index(j) >= primal_size) {
131 return absl::InvalidArgumentError(absl::StrCat(
132 "Variable index of ", i, "th constraint's ", j, "th nonzero is ",
133 con.var_index(j), " which is not in the allowed range [0, ",
134 primal_size, ")"));
135 }
136 nonzeros_by_column[con.var_index(j)]++;
137 }
138 qp.constraint_lower_bounds[i] = con.lower_bound();
139 qp.constraint_upper_bounds[i] = con.upper_bound();
140 if (include_names) {
141 (*qp.constraint_names)[i] = con.name();
142 }
143 }
144 // To reduce peak RAM usage we construct the constraint matrix in-place.
145 // According to the documentation of `SparseMatrix::insert()` it's efficient
146 // to construct a matrix with insert()s as long as reserve() is called first
147 // and the non-zeros are inserted in increasing order of inner index.
148 // The non-zeros in each input constraint may not be sorted so this is only
149 // efficient with column major format.
150 static_assert(qp.constraint_matrix.IsRowMajor == 0, "See comment.");
151 qp.constraint_matrix.reserve(nonzeros_by_column);
152 for (int i = 0; i < dual_size; ++i) {
153 const auto& con = proto.constraint(i);
154 CHECK_EQ(con.var_index_size(), con.coefficient_size())
155 << " in " << i << "th constraint";
156 if (con.var_index_size() != con.coefficient_size()) {
157 return absl::InvalidArgumentError(
158 absl::StrCat(i, "th constraint has ", con.coefficient_size(),
159 " coefficients, expected ", con.var_index_size()));
160 }
161
162 for (int j = 0; j < con.var_index_size(); ++j) {
163 qp.constraint_matrix.insert(i, con.var_index(j)) = con.coefficient(j);
164 }
165 }
166 if (qp.constraint_matrix.outerSize() > 0) {
167 qp.constraint_matrix.makeCompressed();
168 }
169 // We use triplets-based initialization for the objective matrix because the
170 // objective non-zeros may be in arbitrary order in the input.
171 std::vector<Eigen::Triplet<double, int64_t>> triplets;
172 const auto& quadratic = proto.quadratic_objective();
173 if (quadratic.qvar1_index_size() != quadratic.qvar2_index_size() ||
174 quadratic.qvar1_index_size() != quadratic.coefficient_size()) {
175 return absl::InvalidArgumentError(absl::StrCat(
176 "The quadratic objective has ", quadratic.qvar1_index_size(),
177 " qvar1_indices, ", quadratic.qvar2_index_size(),
178 " qvar2_indices, and ", quadratic.coefficient_size(),
179 " coefficients, expected equal numbers."));
180 }
181 if (quadratic.qvar1_index_size() > 0) {
182 qp.objective_matrix.emplace();
183 qp.objective_matrix->setZero(primal_size);
184 }
185
186 for (int i = 0; i < quadratic.qvar1_index_size(); ++i) {
187 const int index1 = quadratic.qvar1_index(i);
188 const int index2 = quadratic.qvar2_index(i);
189 if (index1 < 0 || index2 < 0 || index1 >= primal_size ||
190 index2 >= primal_size) {
191 return absl::InvalidArgumentError(absl::StrCat(
192 "The quadratic objective's ", i, "th nonzero has indices ", index1,
193 " and ", index2, ", which are not both in the expected range [0, ",
194 primal_size, ")"));
195 }
196 if (index1 != index2) {
197 return absl::InvalidArgumentError(absl::StrCat(
198 "The quadratic objective's ", i,
199 "th nonzero has off-diagonal element at (", index1, ", ", index2,
200 "). Only diagonal objective matrices are supported."));
201 }
202 // Note: `QuadraticProgram` has an implicit "1/2" in front of the quadratic
203 // term.
204 qp.objective_matrix->diagonal()[index1] = 2 * quadratic.coefficient(i);
205 }
206 qp.objective_offset = proto.objective_offset();
207 if (proto.maximize()) {
208 qp.objective_offset *= -1;
209 qp.objective_vector *= -1;
210 if (qp.objective_matrix.has_value()) {
211 qp.objective_matrix->diagonal() *= -1;
212 }
214 }
215 return qp;
216}
217
218absl::Status CanFitInMpModelProto(const QuadraticProgram& qp) {
220 qp, std::numeric_limits<int32_t>::max());
221}
222
223namespace internal {
225 const int64_t largest_ok_size) {
226 const int64_t primal_size = qp.variable_lower_bounds.size();
227 const int64_t dual_size = qp.constraint_lower_bounds.size();
228 bool primal_too_big = primal_size > largest_ok_size;
229 if (primal_too_big) {
230 return absl::InvalidArgumentError(absl::StrCat(
231 "Too many variables (", primal_size, ") to index with an int32_t."));
232 }
233 bool dual_too_big = dual_size > largest_ok_size;
234 if (dual_too_big) {
235 return absl::InvalidArgumentError(absl::StrCat(
236 "Too many constraints (", dual_size, ") to index with an int32_t."));
237 }
238 return absl::OkStatus();
239}
240} // namespace internal
241
242absl::StatusOr<MPModelProto> QpToMpModelProto(const QuadraticProgram& qp) {
244 if (qp.objective_scaling_factor == 0) {
245 return absl::InvalidArgumentError(
246 "objective_scaling_factor cannot be zero.");
247 }
248 const int64_t primal_size = qp.variable_lower_bounds.size();
249 const int64_t dual_size = qp.constraint_lower_bounds.size();
250 MPModelProto proto;
251 if (qp.problem_name.has_value() && !qp.problem_name->empty()) {
252 proto.set_name(*qp.problem_name);
253 }
254 proto.set_objective_offset(qp.objective_scaling_factor * qp.objective_offset);
255 if (qp.objective_scaling_factor < 0) {
256 proto.set_maximize(true);
257 } else {
258 proto.set_maximize(false);
259 }
260
261 proto.mutable_variable()->Reserve(primal_size);
262 for (int64_t i = 0; i < primal_size; ++i) {
263 auto* var = proto.add_variable();
264 var->set_lower_bound(qp.variable_lower_bounds[i]);
265 var->set_upper_bound(qp.variable_upper_bounds[i]);
266 var->set_objective_coefficient(qp.objective_scaling_factor *
267 qp.objective_vector[i]);
268 if (qp.variable_names.has_value() && i < qp.variable_names->size()) {
269 const std::string& name = (*qp.variable_names)[i];
270 if (!name.empty()) {
271 var->set_name(name);
272 }
273 }
274 }
275
276 proto.mutable_constraint()->Reserve(dual_size);
277 for (int64_t i = 0; i < dual_size; ++i) {
278 auto* con = proto.add_constraint();
279 con->set_lower_bound(qp.constraint_lower_bounds[i]);
280 con->set_upper_bound(qp.constraint_upper_bounds[i]);
281 if (qp.constraint_names.has_value() && i < qp.constraint_names->size()) {
282 const std::string& name = (*qp.constraint_names)[i];
283 if (!name.empty()) {
284 con->set_name(name);
285 }
286 }
287 }
288
289 using InnerIterator =
290 ::Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>::InnerIterator;
291 for (int64_t col = 0; col < qp.constraint_matrix.cols(); ++col) {
292 for (InnerIterator iter(qp.constraint_matrix, col); iter; ++iter) {
293 auto* con = proto.mutable_constraint(iter.row());
294 // To avoid reallocs during the inserts, we could count the nonzeros
295 // and `reserve()` before filling.
296 con->add_var_index(iter.col());
297 con->add_coefficient(iter.value());
298 }
299 }
300
301 // Some OR tools decide the objective is quadratic based on
302 // `has_quadratic_objective()` rather than on
303 // `quadratic_objective_size() == 0`, so don't create the quadratic objective
304 // for linear programs.
305 if (!IsLinearProgram(qp)) {
306 auto* quadratic_objective = proto.mutable_quadratic_objective();
307 const auto& diagonal = qp.objective_matrix->diagonal();
308 for (int64_t i = 0; i < diagonal.size(); ++i) {
309 if (diagonal[i] != 0.0) {
310 quadratic_objective->add_qvar1_index(i);
311 quadratic_objective->add_qvar2_index(i);
312 // Undo the implicit (1/2) term in `QuadraticProgram`'s objective.
313 quadratic_objective->add_coefficient(qp.objective_scaling_factor *
314 diagonal[i] / 2.0);
315 }
316 }
317 }
318
319 return proto;
320}
321
322std::string ToString(const QuadraticProgram& qp, int64_t max_size) {
323 auto object_name = [](int64_t index,
324 const std::optional<std::vector<std::string>>& names,
325 absl::string_view prefix) {
326 if (names.has_value()) {
327 CHECK_LT(index, names->size());
328 return (*names)[index];
329 }
330 return absl::StrCat(prefix, index);
331 };
332 auto variable_name = [&qp, &object_name](int64_t var_index) {
333 return object_name(var_index, qp.variable_names, "x");
334 };
335 auto constraint_name = [&qp, &object_name](int64_t con_index) {
336 return object_name(con_index, qp.constraint_names, "c");
337 };
338
339 if (auto status = ValidateQuadraticProgramDimensions(qp); !status.ok()) {
340 return absl::StrCat("Quadratic program with inconsistent dimensions: ",
341 status.message());
342 }
343
344 std::string result;
345 if (qp.problem_name.has_value()) {
346 absl::StrAppend(&result, *qp.problem_name, ":\n");
347 }
348 absl::StrAppend(
349 &result, (qp.objective_scaling_factor < 0.0 ? "maximize " : "minimize "),
351 for (int64_t i = 0; i < qp.objective_vector.size(); ++i) {
352 if (qp.objective_vector[i] != 0.0) {
353 absl::StrAppend(&result, " + ", qp.objective_vector[i], " ",
354 variable_name(i));
355 if (result.size() >= max_size) break;
356 }
357 }
358 if (qp.objective_matrix.has_value()) {
359 result.append(" + 1/2 * (");
360 auto obj_diagonal = qp.objective_matrix->diagonal();
361 for (int64_t i = 0; i < obj_diagonal.size(); ++i) {
362 if (obj_diagonal[i] != 0.0) {
363 absl::StrAppend(&result, " + ", obj_diagonal[i], " ", variable_name(i),
364 "^2");
365 if (result.size() >= max_size) break;
366 }
367 }
368 // Closes the objective matrix expression.
369 result.append(")");
370 }
371 // Closes the objective scaling factor expression.
372 result.append(")\n");
373
374 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>
375 constraint_matrix_transpose = qp.constraint_matrix.transpose();
376 for (int64_t constraint_idx = 0;
377 constraint_idx < constraint_matrix_transpose.outerSize();
378 ++constraint_idx) {
379 absl::StrAppend(&result, constraint_name(constraint_idx), ":");
380 if (qp.constraint_lower_bounds[constraint_idx] !=
381 -std::numeric_limits<double>::infinity()) {
382 absl::StrAppend(&result, " ", qp.constraint_lower_bounds[constraint_idx],
383 " <=");
384 }
385 for (decltype(constraint_matrix_transpose)::InnerIterator it(
386 constraint_matrix_transpose, constraint_idx);
387 it; ++it) {
388 absl::StrAppend(&result, " + ", it.value(), " ",
389 variable_name(it.index()));
390 if (result.size() >= max_size) break;
391 }
392 if (qp.constraint_upper_bounds[constraint_idx] !=
393 std::numeric_limits<double>::infinity()) {
394 absl::StrAppend(&result,
395 " <= ", qp.constraint_upper_bounds[constraint_idx]);
396 }
397 result.append("\n");
398 }
399 result.append("Bounds\n");
400 for (int64_t i = 0; i < qp.variable_lower_bounds.size(); ++i) {
401 if (qp.variable_lower_bounds[i] ==
402 -std::numeric_limits<double>::infinity()) {
403 if (qp.variable_upper_bounds[i] ==
404 std::numeric_limits<double>::infinity()) {
405 absl::StrAppend(&result, variable_name(i), " free\n");
406 } else {
407 absl::StrAppend(&result, variable_name(i),
408 " <= ", qp.variable_upper_bounds[i], "\n");
409 }
410 } else {
411 if (qp.variable_upper_bounds[i] ==
412 std::numeric_limits<double>::infinity()) {
413 absl::StrAppend(&result, variable_name(i),
414 " >= ", qp.variable_lower_bounds[i], "\n");
415
416 } else {
417 absl::StrAppend(&result, qp.variable_lower_bounds[i],
418 " <= ", variable_name(i),
419 " <= ", qp.variable_upper_bounds[i], "\n");
420 }
421 }
422 if (result.size() >= max_size) break;
423 }
424 if (result.size() > max_size) {
425 result.resize(max_size - 4);
426 result.append("...\n");
427 }
428 return result;
429}
430
432 std::vector<Eigen::Triplet<double, int64_t>> triplets,
433 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix) {
434 using Triplet = Eigen::Triplet<double, int64_t>;
435 std::sort(triplets.begin(), triplets.end(),
436 [](const Triplet& lhs, const Triplet& rhs) {
437 return std::tie(lhs.col(), lhs.row()) <
438 std::tie(rhs.col(), rhs.row());
439 });
440
441 // The triplets are allowed to contain duplicate entries (and intentionally
442 // do for the diagonals of the objective matrix). For efficiency of insert and
443 // reserve, merge the duplicates first.
445
446 std::vector<int64_t> num_column_entries(matrix.cols());
447 for (const Triplet& triplet : triplets) {
448 ++num_column_entries[triplet.col()];
449 }
450 // NOTE: `reserve()` takes column counts because matrix is in column major
451 // order.
452 matrix.reserve(num_column_entries);
453 for (const Triplet& triplet : triplets) {
454 matrix.insert(triplet.row(), triplet.col()) = triplet.value();
455 }
456 if (matrix.outerSize() > 0) {
457 matrix.makeCompressed();
458 }
459}
460
461namespace internal {
463 std::vector<Eigen::Triplet<double, int64_t>>& triplets) {
464 if (triplets.empty()) return;
465 auto output_iter = triplets.begin();
466 for (auto p = output_iter + 1; p != triplets.end(); ++p) {
467 if (output_iter->row() == p->row() && output_iter->col() == p->col()) {
468 *output_iter = {output_iter->row(), output_iter->col(),
469 output_iter->value() + p->value()};
470 } else {
471 ++output_iter;
472 if (output_iter != p) { // Small optimization - skip no-op copies.
473 *output_iter = *p;
474 }
475 }
476 }
477 // `*output_iter` is the last output value, so erase everything after that.
478 triplets.erase(output_iter + 1, triplets.end());
479}
480} // namespace internal
481} // namespace operations_research::pdlp
IntegerValue size
#define RETURN_IF_ERROR(expr)
CpModelProto proto
The output proto.
const std::string name
A name for logging purposes.
IntVar * var
absl::Status status
Definition g_gurobi.cc:44
int index
ColIndex col
Definition markowitz.cc:187
absl::Status TestableCanFitInMpModelProto(const QuadraticProgram &qp, const int64_t largest_ok_size)
void CombineRepeatedTripletsInPlace(std::vector< Eigen::Triplet< double, int64_t > > &triplets)
Validation utilities for solvers.proto.
absl::StatusOr< MPModelProto > QpToMpModelProto(const QuadraticProgram &qp)
absl::Status ValidateQuadraticProgramDimensions(const QuadraticProgram &qp)
absl::StatusOr< QuadraticProgram > QpFromMpModelProto(const MPModelProto &proto, bool relax_integer_variables, bool include_names)
bool IsLinearProgram(const QuadraticProgram &qp)
std::string ToString(const QuadraticProgram &qp, int64_t max_size)
void SetEigenMatrixFromTriplets(std::vector< Eigen::Triplet< double, int64_t > > triplets, Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix)
absl::Status CanFitInMpModelProto(const QuadraticProgram &qp)
int var_index
Definition search.cc:3268
std::optional< Eigen::DiagonalMatrix< double, Eigen::Dynamic > > objective_matrix
Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > constraint_matrix
std::optional< std::vector< std::string > > constraint_names
std::optional< std::vector< std::string > > variable_names