Google OR-Tools v9.11
a fast and portable software suite for combinatorial optimization
Loading...
Searching...
No Matches
quadratic_program.h
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
14#ifndef PDLP_QUADRATIC_PROGRAM_H_
15#define PDLP_QUADRATIC_PROGRAM_H_
16
17#include <cstdint>
18#include <limits>
19#include <optional>
20#include <string>
21#include <utility>
22#include <vector>
23
24#include "Eigen/Core"
25#include "Eigen/SparseCore"
26#include "absl/status/status.h"
27#include "absl/status/statusor.h"
28#include "ortools/linear_solver/linear_solver.pb.h"
29
31
32// Represents the quadratic program (QP):
33// min_x (`objective_vector`^T x + (1/2) x^T `objective_matrix` x) s.t.
34// `constraint_lower_bounds` <= `constraint_matrix` x
35// <= `constraint_upper_bounds`
36// `variable_lower_bounds` <= x <= `variable_upper_bounds`
37//
38// `constraint_lower_bounds` and `variable_lower_bounds` may include negative
39// infinities. `constraint_upper_bounds` and `variable_upper_bounds` may
40// contain positive infinities. Other than that all entries of all fields must
41// be finite. The `objective_matrix` must be diagonal and non-negative.
42//
43// `variable_lower_bounds`, `variable_upper_bounds`, `objective_vector`,
44// `objective_matrix` (if it has a value), and `variable_names` (if it has a
45// value) must all have the same size as the number of columns in
46// `constraint_matrix`. `constraint_lower_bounds`, `constraint_upper_bounds`,
47// and `constraint_names` (if it has a value) must all have the same size as the
48// number of rows in `constraint_matrix`. Consistency of these values is checked
49// by `ValidateQuadraticProgramDimensions()`.
50//
51// For convenience, the struct also stores `scaling_factor` and
52// `objective_offset`. These factors can be used to transform objective values
53// based on the problem definition above into objective values that are
54// meaningful for the user. See `ApplyObjectiveScalingAndOffset`.
55//
56// This struct is also intended for use with linear programs (LPs), which are
57// QPs with a zero `objective_matrix`.
58//
59// The dual is documented at
60// https://developers.google.com/optimization/lp/pdlp_math.
62 QuadraticProgram(int64_t num_variables, int64_t num_constraints) {
63 ResizeAndInitialize(num_variables, num_constraints);
64 }
66
67 // `QuadraticProgram` may be copied or moved. `Eigen::SparseMatrix` doesn't
68 // have move operations so we use custom implementations based on swap.
69 QuadraticProgram(const QuadraticProgram& other) = default;
71 : objective_vector(std::move(other.objective_vector)),
72 objective_matrix(std::move(other.objective_matrix)),
73 constraint_lower_bounds(std::move(other.constraint_lower_bounds)),
74 constraint_upper_bounds(std::move(other.constraint_upper_bounds)),
75 variable_lower_bounds(std::move(other.variable_lower_bounds)),
76 variable_upper_bounds(std::move(other.variable_upper_bounds)),
77 problem_name(std::move(other.problem_name)),
78 variable_names(std::move(other.variable_names)),
79 constraint_names(std::move(other.constraint_names)),
80 objective_offset(other.objective_offset),
81 objective_scaling_factor(other.objective_scaling_factor) {
82 constraint_matrix.swap(other.constraint_matrix);
83 }
86 objective_vector = std::move(other.objective_vector);
87 objective_matrix = std::move(other.objective_matrix);
88 constraint_matrix.swap(other.constraint_matrix);
89 constraint_lower_bounds = std::move(other.constraint_lower_bounds);
90 constraint_upper_bounds = std::move(other.constraint_upper_bounds);
91 variable_lower_bounds = std::move(other.variable_lower_bounds);
92 variable_upper_bounds = std::move(other.variable_upper_bounds);
93 problem_name = std::move(other.problem_name);
94 variable_names = std::move(other.variable_names);
95 constraint_names = std::move(other.constraint_names);
96 objective_offset = other.objective_offset;
97 objective_scaling_factor = other.objective_scaling_factor;
98 return *this;
99 }
100
101 // Initializes the quadratic program with `num_variables` variables and
102 // `num_constraints` constraints. Lower and upper bounds are set to negative
103 // and positive infinity, respectively. `objective_matrix` is cleared. All
104 // other matrices and vectors are set to zero. Resets the optional names
105 // (`program_name`, `variable_names`, and `constraint_names`).
106 // `objective_offset` is set to 0 and `objective_scaling_factor` is set to 1.
107 void ResizeAndInitialize(int64_t num_variables, int64_t num_constraints) {
108 constexpr double kInfinity = std::numeric_limits<double>::infinity();
109 objective_vector = Eigen::VectorXd::Zero(num_variables);
110 objective_matrix.reset();
111 constraint_matrix.resize(num_constraints, num_variables);
113 Eigen::VectorXd::Constant(num_constraints, -kInfinity);
115 Eigen::VectorXd::Constant(num_constraints, kInfinity);
117 Eigen::VectorXd::Constant(num_variables, -kInfinity);
118 variable_upper_bounds = Eigen::VectorXd::Constant(num_variables, kInfinity);
119 problem_name.reset();
120 variable_names.reset();
121 constraint_names.reset();
122 objective_offset = 0.0;
124 }
125
126 // Returns `objective_scaling_factor * (objective + objective_offset)`.
127 // `objective_scaling_factor` is useful for modeling maximization problems.
128 // For example, max c^T x = -1 * min (-c)^T x. `objective_offset` can be a
129 // by-product of presolve transformations that eliminate variables.
130 double ApplyObjectiveScalingAndOffset(double objective) const {
131 return objective_scaling_factor * (objective + objective_offset);
132 }
133
134 Eigen::VectorXd objective_vector;
135 // If this field isn't set, the `objective matrix` is interpreted to be zero,
136 // i.e., this is a linear programming problem.
137 std::optional<Eigen::DiagonalMatrix<double, Eigen::Dynamic>> objective_matrix;
138 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> constraint_matrix;
141
142 std::optional<std::string> problem_name;
143 std::optional<std::vector<std::string>> variable_names;
144 std::optional<std::vector<std::string>> constraint_names;
145
146 // These fields are provided for convenience; they don't change the
147 // mathematical definition of the problem, but they change the objective
148 // values reported to the user.
151};
152
153// Returns `InvalidArgumentError` if vector or matrix dimensions are
154// inconsistent. Returns `OkStatus` otherwise.
156
157inline bool IsLinearProgram(const QuadraticProgram& qp) {
158 return !qp.objective_matrix.has_value();
159}
160
161// Converts an `MPModelProto` into a `QuadraticProgram`.
162// Returns an error if general constraints are present.
163// If `relax_integer_variables` is true integer variables are relaxed to
164// continuous; otherwise integer variables are an error.
165// If `include_names` is true (the default is false), the problem, constraint,
166// and variable names are included in the `QuadraticProgram`; otherwise they are
167// left empty.
168// Maximization problems are converted to minimization by negating the
169// objective and setting `objective_scaling_factor` to -1, which preserves the
170// reported objective values.
171absl::StatusOr<QuadraticProgram> QpFromMpModelProto(
172 const MPModelProto& proto, bool relax_integer_variables,
173 bool include_names = false);
174
175// Returns `InvalidArgumentError` if `qp` is too large to convert to
176// `MPModelProto` and `OkStatus` otherwise.
177absl::Status CanFitInMpModelProto(const QuadraticProgram& qp);
178
179// Converts a `QuadraticProgram` into an `MPModelProto`. To preserve objective
180// values in the conversion, `objective_vector`, `objective_matrix`, and
181// `objective_offset` are scaled by `objective_scaling_factor`, and if
182// `objective_scaling_factor` is negative, then the proto is a maximization
183// problem (otherwise it's a minimization problem). Returns
184// `InvalidArgumentError` if `objective_scaling_factor` is zero or if
185// `CanFitInMpModelProto()` fails.
186absl::StatusOr<MPModelProto> QpToMpModelProto(const QuadraticProgram& qp);
187
188// Returns a "pretty" version of `qp`, truncating to at most `max_size`
189// characters. This is for debugging purposes only - the format may change
190// without notice. Although this output is vaguely similar to "LP format", it is
191// not actually compatible with "LP format".
192std::string ToString(const QuadraticProgram& qp, int64_t max_size = 1'000'000);
193
194// Like `matrix.setFromTriplets(triplets)`, except that `setFromTriplets`
195// results in having three copies of the nonzeros in memory at the same time,
196// because it first fills one matrix from triplets, and then transposes it into
197// another. This avoids having the third copy in memory by sorting the triplets,
198// reserving space in the matrix, and then inserting in sorted order.
199// Compresses the matrix (`SparseMatrix.makeCompressed()`) after loading it.
200// NOTE: This intentionally passes `triplets` by value, because it modifies
201// them. To avoid the copy, pass a move reference.
203 std::vector<Eigen::Triplet<double, int64_t>> triplets,
204 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix);
205
206// Utility functions for internal use only.
207namespace internal {
208// Like `CanFitInMpModelProto()` but has an extra argument for the largest
209// number of variables, constraints, or objective non-zeros that should be
210// counted as convertible. `CanFitInMpModelProto()` passes 2^31 - 1 for this
211// argument and unit tests pass small values.
212absl::Status TestableCanFitInMpModelProto(const QuadraticProgram& qp,
213 int64_t largest_ok_size);
214
215// Modifies `triplets` in place, combining consecutive entries with the same row
216// and column, summing their values. This is most effective if `triplets` are
217// sorted by row and column, so that multiple entries for the same entry will be
218// consecutive.
220 std::vector<Eigen::Triplet<double, int64_t>>& triplets);
221} // namespace internal
222} // namespace operations_research::pdlp
223
224#endif // PDLP_QUADRATIC_PROGRAM_H_
CpModelProto proto
The output proto.
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)
double ApplyObjectiveScalingAndOffset(double objective) const
QuadraticProgram & operator=(QuadraticProgram &&other)
QuadraticProgram(const QuadraticProgram &other)=default
std::optional< Eigen::DiagonalMatrix< double, Eigen::Dynamic > > objective_matrix
Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > constraint_matrix
QuadraticProgram & operator=(const QuadraticProgram &other)=default
std::optional< std::vector< std::string > > constraint_names
QuadraticProgram(int64_t num_variables, int64_t num_constraints)
std::optional< std::vector< std::string > > variable_names
void ResizeAndInitialize(int64_t num_variables, int64_t num_constraints)
QuadraticProgram(QuadraticProgram &&other) noexcept