Google OR-Tools v9.11
a fast and portable software suite for combinatorial optimization
|
#include <quadratic_program.h>
Public Member Functions | |
QuadraticProgram (int64_t num_variables, int64_t num_constraints) | |
QuadraticProgram () | |
QuadraticProgram (const QuadraticProgram &other)=default | |
QuadraticProgram (QuadraticProgram &&other) noexcept | |
QuadraticProgram & | operator= (const QuadraticProgram &other)=default |
QuadraticProgram & | operator= (QuadraticProgram &&other) |
void | ResizeAndInitialize (int64_t num_variables, int64_t num_constraints) |
double | ApplyObjectiveScalingAndOffset (double objective) const |
Public Attributes | |
Eigen::VectorXd | objective_vector |
std::optional< Eigen::DiagonalMatrix< double, Eigen::Dynamic > > | objective_matrix |
Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > | constraint_matrix |
Eigen::VectorXd | constraint_lower_bounds |
Eigen::VectorXd | constraint_upper_bounds |
Eigen::VectorXd | variable_lower_bounds |
Eigen::VectorXd | variable_upper_bounds |
std::optional< std::string > | problem_name |
std::optional< std::vector< std::string > > | variable_names |
std::optional< std::vector< std::string > > | constraint_names |
double | objective_offset |
double | objective_scaling_factor |
Represents the quadratic program (QP): min_x (objective_vector
^T x + (1/2) x^T objective_matrix
x) s.t. constraint_lower_bounds
<= constraint_matrix
x <= constraint_upper_bounds
variable_lower_bounds
<= x <= variable_upper_bounds
constraint_lower_bounds
and variable_lower_bounds
may include negative infinities. constraint_upper_bounds
and variable_upper_bounds
may contain positive infinities. Other than that all entries of all fields must be finite. The objective_matrix
must be diagonal and non-negative.
variable_lower_bounds
, variable_upper_bounds
, objective_vector
, objective_matrix
(if it has a value), and variable_names
(if it has a value) must all have the same size as the number of columns in constraint_matrix
. constraint_lower_bounds
, constraint_upper_bounds
, and constraint_names
(if it has a value) must all have the same size as the number of rows in constraint_matrix
. Consistency of these values is checked by ValidateQuadraticProgramDimensions()
.
For convenience, the struct also stores scaling_factor
and objective_offset
. These factors can be used to transform objective values based on the problem definition above into objective values that are meaningful for the user. See ApplyObjectiveScalingAndOffset
.
This struct is also intended for use with linear programs (LPs), which are QPs with a zero objective_matrix
.
The dual is documented at https://developers.google.com/optimization/lp/pdlp_math.
Definition at line 61 of file quadratic_program.h.
|
inline |
Definition at line 62 of file quadratic_program.h.
|
inline |
Definition at line 65 of file quadratic_program.h.
|
default |
QuadraticProgram
may be copied or moved. Eigen::SparseMatrix
doesn't have move operations so we use custom implementations based on swap.
|
inlinenoexcept |
Definition at line 70 of file quadratic_program.h.
|
inline |
Returns objective_scaling_factor * (objective + objective_offset)
. objective_scaling_factor
is useful for modeling maximization problems. For example, max c^T x = -1 * min (-c)^T x. objective_offset
can be a by-product of presolve transformations that eliminate variables.
Definition at line 130 of file quadratic_program.h.
|
default |
|
inline |
Definition at line 85 of file quadratic_program.h.
|
inline |
Initializes the quadratic program with num_variables
variables and num_constraints
constraints. Lower and upper bounds are set to negative and positive infinity, respectively. objective_matrix
is cleared. All other matrices and vectors are set to zero. Resets the optional names (program_name
, variable_names
, and constraint_names
). objective_offset
is set to 0 and objective_scaling_factor
is set to 1.
Definition at line 107 of file quadratic_program.h.
Eigen::VectorXd operations_research::pdlp::QuadraticProgram::constraint_lower_bounds |
Definition at line 139 of file quadratic_program.h.
Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> operations_research::pdlp::QuadraticProgram::constraint_matrix |
Definition at line 138 of file quadratic_program.h.
std::optional<std::vector<std::string> > operations_research::pdlp::QuadraticProgram::constraint_names |
Definition at line 144 of file quadratic_program.h.
Eigen::VectorXd operations_research::pdlp::QuadraticProgram::constraint_upper_bounds |
Definition at line 139 of file quadratic_program.h.
std::optional<Eigen::DiagonalMatrix<double, Eigen::Dynamic> > operations_research::pdlp::QuadraticProgram::objective_matrix |
If this field isn't set, the objective matrix
is interpreted to be zero, i.e., this is a linear programming problem.
Definition at line 137 of file quadratic_program.h.
double operations_research::pdlp::QuadraticProgram::objective_offset |
These fields are provided for convenience; they don't change the mathematical definition of the problem, but they change the objective values reported to the user.
Definition at line 149 of file quadratic_program.h.
double operations_research::pdlp::QuadraticProgram::objective_scaling_factor |
Definition at line 150 of file quadratic_program.h.
Eigen::VectorXd operations_research::pdlp::QuadraticProgram::objective_vector |
Definition at line 134 of file quadratic_program.h.
std::optional<std::string> operations_research::pdlp::QuadraticProgram::problem_name |
Definition at line 142 of file quadratic_program.h.
Eigen::VectorXd operations_research::pdlp::QuadraticProgram::variable_lower_bounds |
Definition at line 140 of file quadratic_program.h.
std::optional<std::vector<std::string> > operations_research::pdlp::QuadraticProgram::variable_names |
Definition at line 143 of file quadratic_program.h.
Eigen::VectorXd operations_research::pdlp::QuadraticProgram::variable_upper_bounds |
Definition at line 140 of file quadratic_program.h.