OR-Tools  9.3
sharded_optimization_utils.h
Go to the documentation of this file.
1// Copyright 2010-2021 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// These are internal helper functions and classes that implicitly or explicitly
15// operate on a ShardedQuadraticProgram. Utilities that are purely linear
16// algebra operations (e.g., norms) should be defined in sharder.h instead.
17
18#ifndef PDLP_SHARDED_OPTIMIZATION_UTILS_H_
19#define PDLP_SHARDED_OPTIMIZATION_UTILS_H_
20
21#include <limits>
22#include <random>
23
24#include "Eigen/Core"
25#include "absl/types/optional.h"
28#include "ortools/pdlp/solve_log.pb.h"
29
31
32// This computes weighted averages of vectors.
33// It satisfies the following: if all the averaged vectors have component i
34// equal to x then the average has component i exactly equal to x, without any
35// floating-point roundoff. In practice the above is probably still true with
36// "equal to x" replaced with "at least x" or "at most x". However unrealistic
37// counter examples probably exist involving a new item with weight 10^15 times
38// greater than the total weight so far.
40 public:
41 // Initializes the weighted average by creating a vector sized according to
42 // the number of elements in the sharder. Retains the pointer to sharder, so
43 // the sharder must outlive this object.
44 explicit ShardedWeightedAverage(const Sharder* sharder);
45
48
49 // Adds the datapoint to the average with the given weight. CHECK-fails if
50 // the weight is negative.
51 void Add(const Eigen::VectorXd& datapoint, double weight);
52
53 // Clears the sum to zero, i.e., just constructed.
54 void Clear();
55
56 // Returns true if there is at least one term in the average with a positive
57 // weight.
58 bool HasNonzeroWeight() const { return sum_weights_ > 0.0; }
59
60 // Returns the sum of the weights of the datapoints added so far.
61 double Weight() const { return sum_weights_; }
62
63 // Computes the weighted average of the datapoints added so far, i.e.,
64 // sum_i weight[i] * datapoint[i] / sum_i weight[i]. The results are set to
65 // zero if HasNonzeroWeight() is false.
66 Eigen::VectorXd ComputeAverage() const;
67
68 int NumTerms() const { return num_terms_; }
69
70 private:
71 Eigen::VectorXd average_;
72 double sum_weights_ = 0.0;
73 int num_terms_ = 0;
74 const Sharder* sharder_;
75};
76
77// Returns a QuadraticProgramStats object for a ShardedQuadraticProgram.
78QuadraticProgramStats ComputeStats(const ShardedQuadraticProgram& qp,
79 double infinite_constraint_bound_threshold =
80 std::numeric_limits<double>::infinity());
81
82// LInfRuizRescaling and L2NormRescaling rescale the (scaled) constraint matrix
83// of the LP by updating the scaling vectors in-place. More specifically, the
84// (scaled) constraint matrix always has the format: diag(row_scaling_vec) *
85// sharded_qp.Qp().constraint_matrix * diag(col_scaling_vec), and
86// row_scaling_vec and col_scaling_vec are updated in-place. On input,
87// row_scaling_vec and col_scaling_vec provide the initial scaling.
88
89// With each iteration of LInfRuizRescaling scaling, row_scaling_vec
90// (col_scaling_vec) is divided by the sqrt of the row (col) LInf norm of the
91// current (scaled) constraint matrix. The (scaled) constraint matrix approaches
92// having all row and column LInf norms equal to 1 as the number of iterations
93// goes to infinity. This convergence is fast (linear). More details of Ruiz
94// rescaling algorithm can be found at:
95// http://www.numerical.rl.ac.uk/reports/drRAL2001034.pdf.
97 const int num_iterations,
98 Eigen::VectorXd& row_scaling_vec,
99 Eigen::VectorXd& col_scaling_vec);
100
101// L2NormRescaling divides row_scaling_vec (col_scaling_vec) by the sqrt of the
102// row (col) L2 norm of the current (scaled) constraint matrix. Unlike
103// LInfRescaling, this function does only one iteration because the scaling
104// procedure does not converge in general. This is not Ruiz rescaling for the L2
105// norm.
107 Eigen::VectorXd& row_scaling_vec,
108 Eigen::VectorXd& col_scaling_vec);
109
113};
114
116 Eigen::VectorXd row_scaling_vec;
117 Eigen::VectorXd col_scaling_vec;
118};
119
120// Applies the rescaling specified by rescaling_options to sharded_qp (in
121// place). Returns the scaling vectors that were applied.
122ScalingVectors ApplyRescaling(const RescalingOptions& rescaling_options,
123 ShardedQuadraticProgram& sharded_qp);
124
126 double value = 0.0;
127 Eigen::VectorXd gradient;
128};
129
130// Computes the value of the primal part of the Lagrangian function defined at
131// https://developers.google.com/optimization/lp/pdlp_math, i.e.,
132// c'x + (1/2) x'Qx - y'Ax and its gradient with respect to the primal variables
133// x, i.e., c + Qx - A'y. The dual_product argument is A'y. Note: The objective
134// constant is omitted. The result is undefined and invalid if any primal bounds
135// are violated.
137 const Eigen::VectorXd& primal_solution,
138 const Eigen::VectorXd& dual_product);
139
140// Returns a subderivative of the concave dual penalty function that appears in
141// the Lagrangian: -p(dual; -constraint_upper_bound, -constraint_lower_bound) =
142// { constraint_upper_bound * dual when dual < 0, 0 when dual == 0, and
143// constraint_lower_bound * dual when dual > 0}
144// (as defined at https://developers.google.com/optimization/lp/pdlp_math).
145// The subderivative is not necessarily unique when dual == 0. In this case, if
146// only one of the bounds is finite, we return that one. If both are finite, we
147// return the primal product projected onto the bounds, which causes the dual
148// Lagrangian gradient to be zero when the constraint is not violated. If both
149// are infinite, we return zero. The value returned is valid only when the
150// function is finite-valued.
151double DualSubgradientCoefficient(const double constraint_lower_bound,
152 const double constraint_upper_bound,
153 const double dual,
154 const double primal_product);
155
156// Computes the value of the dual part of the Lagrangian function defined at
157// https://developers.google.com/optimization/lp/pdlp_math, i.e., -h^*(y) and
158// the gradient of the Lagrangian with respect to the dual variables y, i.e.,
159// -Ax - \grad_y h^*(y). Note the asymmetry with ComputePrimalGradient: the term
160// -y'Ax is not part of the value. Because h^*(y) is piece-wise linear, a
161// subgradient is returned at a point of non- smoothness. The primal_product
162// argument is Ax. The result is undefined and invalid if any duals violate
163// their bounds.
165 const Eigen::VectorXd& dual_solution,
166 const Eigen::VectorXd& primal_product);
167
172};
173
174// Estimates the maximum singular value of A by applying the method of power
175// iteration to A^T A. If a primal and/or dual solution is provided, restricts
176// to the "active" part of A, that is, the columns (rows) for variables that are
177// not at their bounds in the solution. The estimate will have
178// desired_relative_error with probability at least 1 - failure_probability.
179// The number of iterations will be approximately
180// log(primal_size / failure_probability^2)/(2 * desired_relative_error).
181// Uses a mersenne-twister portable random number generator to generate the
182// starting point for the power method, in order to have deterministic results.
184 const ShardedQuadraticProgram& sharded_qp,
185 const absl::optional<Eigen::VectorXd>& primal_solution,
186 const absl::optional<Eigen::VectorXd>& dual_solution,
187 const double desired_relative_error, const double failure_probability,
188 std::mt19937& mt_generator);
189
190// Checks if the lower and upper bounds of the problem are consistent, i.e. for
191// each variable and constraint bound we have lower_bound <= upper_bound. If
192// the input is consistent the method returns true, otherwise it returns false.
193// See also HasValidBounds(const QuadraticProgram&).
194bool HasValidBounds(const ShardedQuadraticProgram& sharded_qp);
195
196// Projects a primal vector onto the variable bounds constraints.
198 Eigen::VectorXd& primal);
199
200// Projects the dual vector to the dual variable bounds; see
201// https://developers.google.com/optimization/lp/pdlp_math#dual_variable_bounds.
203 Eigen::VectorXd& dual);
204
205} // namespace operations_research::pdlp
206
207#endif // PDLP_SHARDED_OPTIMIZATION_UTILS_H_
void Add(const Eigen::VectorXd &datapoint, double weight)
ShardedWeightedAverage & operator=(ShardedWeightedAverage &&)=default
ShardedWeightedAverage(ShardedWeightedAverage &&)=default
LagrangianPart ComputePrimalGradient(const ShardedQuadraticProgram &sharded_qp, const VectorXd &primal_solution, const VectorXd &dual_product)
void LInfRuizRescaling(const ShardedQuadraticProgram &sharded_qp, const int num_iterations, VectorXd &row_scaling_vec, VectorXd &col_scaling_vec)
LagrangianPart ComputeDualGradient(const ShardedQuadraticProgram &sharded_qp, const Eigen::VectorXd &dual_solution, const Eigen::VectorXd &primal_product)
double DualSubgradientCoefficient(const double constraint_lower_bound, const double constraint_upper_bound, const double dual, const double primal_product)
bool HasValidBounds(const QuadraticProgram &qp)
SingularValueAndIterations EstimateMaximumSingularValueOfConstraintMatrix(const ShardedQuadraticProgram &sharded_qp, const absl::optional< VectorXd > &primal_solution, const absl::optional< VectorXd > &dual_solution, const double desired_relative_error, const double failure_probability, std::mt19937 &mt_generator)
void ProjectToDualVariableBounds(const ShardedQuadraticProgram &sharded_qp, VectorXd &dual)
void L2NormRescaling(const ShardedQuadraticProgram &sharded_qp, VectorXd &row_scaling_vec, VectorXd &col_scaling_vec)
ScalingVectors ApplyRescaling(const RescalingOptions &rescaling_options, ShardedQuadraticProgram &sharded_qp)
QuadraticProgramStats ComputeStats(const ShardedQuadraticProgram &qp, const double infinite_constraint_bound_threshold)
void ProjectToPrimalVariableBounds(const ShardedQuadraticProgram &sharded_qp, VectorXd &primal)
int64_t weight
Definition: pack.cc:510