OR-Tools  9.3
iteration_stats.cc
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
15
16#include <algorithm>
17#include <cmath>
18#include <cstdint>
19#include <limits>
20#include <random>
21#include <string>
22#include <utility>
23#include <vector>
24
25#include "Eigen/Core"
26#include "Eigen/SparseCore"
27#include "absl/random/distributions.h"
28#include "absl/strings/str_cat.h"
29#include "absl/strings/str_format.h"
30#include "absl/strings/string_view.h"
31#include "absl/types/optional.h"
37#include "ortools/pdlp/solve_log.pb.h"
38#include "ortools/pdlp/solvers.pb.h"
39
41namespace {
42
43using ::Eigen::VectorXd;
44
45// ResidualNorms contains four measures of the infeasibility of a primal or
46// dual solution. "objective_correction" is the (additive) adjustment to the
47// objective function from the reduced costs. "objective_full_correction" is the
48// (additive) adjustment to the objective function if all dual residuals were
49// set to zero, while l_inf_residual and l_2_residual are the L_infinity and L_2
50// norms of the residuals (portions of the primal gradient not included in the
51// reduced costs).
52struct ResidualNorms {
57};
58
59// Computes norms of the primal residual infeasibilities (b - A x) of the
60// unscaled problem. Note the primal residuals of the unscaled problem are equal
61// to those of the scaled problem divided by row_scaling_vec. Does not perform
62// any corrections (so the returned .correction == 0.0). sharded_qp is assumed
63// to be the scaled problem. If use_homogeneous_primal_bounds is set to true
64// the residuals are computed with upper and lower bounds zeroed out (note that
65// we only zero out the bounds that are finite in the original problem).
66ResidualNorms PrimalResidualNorms(
67 const ShardedQuadraticProgram& sharded_qp, const VectorXd& row_scaling_vec,
68 const VectorXd& scaled_primal_solution,
69 bool use_homogeneous_constraint_bounds = false) {
70 const QuadraticProgram& qp = sharded_qp.Qp();
71 CHECK_EQ(row_scaling_vec.size(), sharded_qp.DualSize());
72 CHECK_EQ(scaled_primal_solution.size(), sharded_qp.PrimalSize());
73
74 VectorXd primal_product = TransposedMatrixVectorProduct(
75 sharded_qp.TransposedConstraintMatrix(), scaled_primal_solution,
76 sharded_qp.TransposedConstraintMatrixSharder());
77 VectorXd local_l_inf_residual(sharded_qp.DualSharder().NumShards());
78 VectorXd local_sumsq_residual(sharded_qp.DualSharder().NumShards());
79 sharded_qp.DualSharder().ParallelForEachShard(
80 [&](const Sharder::Shard& shard) {
81 const auto lower_bound_shard = shard(qp.constraint_lower_bounds);
82 const auto upper_bound_shard = shard(qp.constraint_upper_bounds);
83 const auto row_scaling_shard = shard(row_scaling_vec);
84 const auto primal_product_shard = shard(primal_product);
85 double l_inf_residual = 0.0;
86 double sumsq_residual = 0.0;
87 for (int64_t i = 0; i < primal_product_shard.size(); ++i) {
88 const double upper_bound = (use_homogeneous_constraint_bounds &&
89 std::isfinite(upper_bound_shard[i]))
90 ? 0.0
91 : upper_bound_shard[i];
92 const double lower_bound = (use_homogeneous_constraint_bounds &&
93 std::isfinite(lower_bound_shard[i]))
94 ? 0.0
95 : lower_bound_shard[i];
96 double scaled_residual = 0.0;
97 if (primal_product_shard[i] > upper_bound) {
98 scaled_residual = primal_product_shard[i] - upper_bound;
99
100 } else if (primal_product_shard[i] < lower_bound) {
101 scaled_residual = lower_bound - primal_product_shard[i];
102 }
103 const double residual = scaled_residual / row_scaling_shard[i];
105 sumsq_residual += residual * residual;
106 }
107 local_l_inf_residual[shard.Index()] = l_inf_residual;
108 local_sumsq_residual[shard.Index()] = sumsq_residual;
109 });
110 return ResidualNorms{
111 .objective_correction = 0.0,
112 .objective_full_correction = 0.0,
113 .l_inf_residual = local_l_inf_residual.lpNorm<Eigen::Infinity>(),
114 .l_2_residual = std::sqrt(local_sumsq_residual.sum()),
115 };
116}
117
118// Decides whether a primal gradient term should be handled as a reduced cost or
119// as a dual residual.
120bool HandlePrimalGradientTermAsReducedCost(double primal_gradient,
121 double primal_value,
122 double lower_bound,
123 double upper_bound) {
124 if (primal_gradient == 0.0) return true;
125 return std::abs(primal_value -
126 (primal_gradient > 0.0 ? lower_bound : upper_bound)) <=
127 std::abs(primal_value);
128}
129
130// Computes norms of the dual residuals and reduced costs of the unscaled
131// problem. Note the primal gradient of the unscaled problem is equal to the
132// scaled primal gradient divided by col_scaling_vec. sharded_qp is assumed to
133// be the scaled problem. See
134// https://developers.google.com/optimization/lp/pdlp_math for details and
135// notation. Primal gradients that have corresponding (finite) bounds (the
136// finite terms from (l^v)^T[r]_+ − (u^v)^T[r]_− in the dual objective), and
137// have |x - b| <= |x| (where x is the variable's value and b is the
138// corresponding bound) are treated as reduced costs and accumulated in
139// objective_correction, while the other primal gradient terms are handled as
140// residual infeasibilities in l_inf_residual and l_2_residual.
141ResidualNorms DualResidualNorms(const ShardedQuadraticProgram& sharded_qp,
142 const VectorXd& col_scaling_vec,
143 const VectorXd& scaled_primal_solution,
144 const VectorXd& scaled_primal_gradient) {
145 const QuadraticProgram& qp = sharded_qp.Qp();
146 CHECK_EQ(col_scaling_vec.size(), sharded_qp.PrimalSize());
147 CHECK_EQ(scaled_primal_gradient.size(), sharded_qp.PrimalSize());
148 VectorXd local_dual_correction(sharded_qp.PrimalSharder().NumShards());
149 VectorXd local_dual_full_correction(sharded_qp.PrimalSharder().NumShards());
150 VectorXd local_l_inf_residual(sharded_qp.PrimalSharder().NumShards());
151 VectorXd local_sumsq_residual(sharded_qp.PrimalSharder().NumShards());
152 sharded_qp.PrimalSharder().ParallelForEachShard(
153 [&](const Sharder::Shard& shard) {
154 const auto lower_bound_shard = shard(qp.variable_lower_bounds);
155 const auto upper_bound_shard = shard(qp.variable_upper_bounds);
156 const auto primal_gradient_shard = shard(scaled_primal_gradient);
157 const auto col_scaling_shard = shard(col_scaling_vec);
158 const auto primal_solution_shard = shard(scaled_primal_solution);
159 double dual_correction = 0.0;
160 double dual_full_correction = 0.0;
161 double l_inf_residual = 0.0;
162 double sumsq_residual = 0.0;
163 for (int64_t i = 0; i < primal_gradient_shard.size(); ++i) {
164 // The corrections use the scaled values because
165 // unscaled_lower_bound = lower_bound * scale and
166 // unscaled_primal_gradient = primal_gradient / scale, so the scales
167 // cancel out.
168 if (primal_gradient_shard[i] == 0.0) continue;
169 const double bound_for_rc = primal_gradient_shard[i] > 0.0
170 ? lower_bound_shard[i]
171 : upper_bound_shard[i];
172 dual_full_correction += bound_for_rc * primal_gradient_shard[i];
173 if (HandlePrimalGradientTermAsReducedCost(
174 primal_gradient_shard[i], primal_solution_shard[i],
175 lower_bound_shard[i], upper_bound_shard[i])) {
176 dual_correction += bound_for_rc * primal_gradient_shard[i];
177 } else {
178 const double scaled_residual = std::abs(primal_gradient_shard[i]);
179 const double residual = scaled_residual / col_scaling_shard[i];
181 sumsq_residual += residual * residual;
182 }
183 }
184 local_dual_correction[shard.Index()] = dual_correction;
185 local_dual_full_correction[shard.Index()] = dual_full_correction;
186 local_l_inf_residual[shard.Index()] = l_inf_residual;
187 local_sumsq_residual[shard.Index()] = sumsq_residual;
188 });
189 return ResidualNorms{
190 .objective_correction = local_dual_correction.sum(),
191 .objective_full_correction = local_dual_full_correction.sum(),
192 .l_inf_residual = local_l_inf_residual.lpNorm<Eigen::Infinity>(),
193 .l_2_residual = std::sqrt(local_sumsq_residual.sum()),
194 };
195}
196
197// Returns Qx.
198VectorXd ObjectiveProduct(const ShardedQuadraticProgram& sharded_qp,
199 const VectorXd& primal_solution) {
200 CHECK_EQ(primal_solution.size(), sharded_qp.PrimalSize());
201 VectorXd result(primal_solution.size());
202 if (IsLinearProgram(sharded_qp.Qp())) {
203 result.setZero();
204 } else {
205 sharded_qp.PrimalSharder().ParallelForEachShard(
206 [&](const Sharder::Shard& shard) {
207 shard(result) =
208 shard(*sharded_qp.Qp().objective_matrix) * shard(primal_solution);
209 });
210 }
211 return result;
212}
213
214// Returns 1/2 x^T Q x (the quadratic term in the objective).
215double QuadraticObjective(const ShardedQuadraticProgram& sharded_qp,
216 const VectorXd& primal_solution,
217 const VectorXd& objective_product) {
218 CHECK_EQ(primal_solution.size(), sharded_qp.PrimalSize());
219 CHECK_EQ(objective_product.size(), sharded_qp.PrimalSize());
220 return 0.5 *
221 Dot(objective_product, primal_solution, sharded_qp.PrimalSharder());
222}
223
224// Returns objective_product + c − A^T y when use_zero_primal_objective =
225// false, and returns − A^T y when use_zero_primal_objective = true.
226// objective_product is passed by copy, and modified in place.
227VectorXd PrimalGradientFromObjectiveProduct(
228 const ShardedQuadraticProgram& sharded_qp, const VectorXd& dual_solution,
229 VectorXd objective_product, bool use_zero_primal_objective = false) {
230 const QuadraticProgram& qp = sharded_qp.Qp();
231 CHECK_EQ(dual_solution.size(), sharded_qp.DualSize());
232 CHECK_EQ(objective_product.size(), sharded_qp.PrimalSize());
233
234 // Note that this modifies objective_product, replacing its entries with
235 // the primal gradient.
236 sharded_qp.ConstraintMatrixSharder().ParallelForEachShard(
237 [&](const Sharder::Shard& shard) {
238 if (use_zero_primal_objective) {
239 shard(objective_product) =
240 -shard(qp.constraint_matrix).transpose() * dual_solution;
241 } else {
242 shard(objective_product) +=
243 shard(qp.objective_vector) -
244 shard(qp.constraint_matrix).transpose() * dual_solution;
245 }
246 });
247 return objective_product;
248}
249
250// Returns the value of y term in the objective of the dual problem, see
251// (l^c)^T[y]_+ − (u^c)^T[y]_− in the dual objective from
252// https://developers.google.com/optimization/lp/pdlp_math.
253double DualObjectiveBoundsTerm(const ShardedQuadraticProgram& sharded_qp,
254 const VectorXd& dual_solution) {
255 const QuadraticProgram& qp = sharded_qp.Qp();
256 return sharded_qp.DualSharder().ParallelSumOverShards(
257 [&](const Sharder::Shard& shard) {
258 // This assumes that the dual variables are feasible, that is, that
259 // the term corresponding to the "y" variables in the dual objective
260 // in https://developers.google.com/optimization/lp/pdlp_math is finite.
261 const auto lower_bound_shard = shard(qp.constraint_lower_bounds);
262 const auto upper_bound_shard = shard(qp.constraint_upper_bounds);
263 const auto dual_shard = shard(dual_solution);
264 // Can't use .dot(.cwiseMin()) because that gives 0 * inf = NaN.
265 double sum = 0.0;
266 for (int64_t i = 0; i < dual_shard.size(); ++i) {
267 if (dual_shard[i] > 0.0) {
268 sum += lower_bound_shard[i] * dual_shard[i];
269 } else if (dual_shard[i] < 0.0) {
270 sum += upper_bound_shard[i] * dual_shard[i];
271 }
272 }
273 return sum;
274 });
275}
276
277// Computes the projection of the vector onto a pseudo-random vector determined
278// by seed_generator. seed_generator is used as the source of a random seed for
279// each shard's portion of the vector.
280double RandomProjection(const VectorXd& vector, const Sharder& sharder,
281 std::mt19937& seed_generator) {
282 std::vector<std::mt19937> shard_seeds;
283 shard_seeds.reserve(sharder.NumShards());
284 for (int shard = 0; shard < sharder.NumShards(); ++shard) {
285 shard_seeds.emplace_back((seed_generator)());
286 }
287 // Computes vector * gaussian_random_vector and
288 // ||gaussian_random_vector||^2 to normalize by afterwards.
289 VectorXd dot_product(sharder.NumShards());
290 VectorXd gaussian_norm_squared(sharder.NumShards());
291 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
292 const auto vector_shard = shard(vector);
293 double shard_dot_product = 0.0;
294 double shard_norm_squared = 0.0;
295 std::mt19937 random{shard_seeds[shard.Index()]};
296 for (int64_t i = 0; i < vector_shard.size(); ++i) {
297 const double projection_element = absl::Gaussian(random, 0.0, 1.0);
298 shard_dot_product += projection_element * vector_shard[i];
299 shard_norm_squared += MathUtil::Square(projection_element);
300 }
301 dot_product[shard.Index()] = shard_dot_product;
302 gaussian_norm_squared[shard.Index()] = shard_norm_squared;
303 });
304 return dot_product.sum() / std::sqrt(gaussian_norm_squared.sum());
305}
306} // namespace
307
308ConvergenceInformation ComputeConvergenceInformation(
309 const ShardedQuadraticProgram& scaled_sharded_qp,
310 const Eigen::VectorXd& col_scaling_vec,
311 const Eigen::VectorXd& row_scaling_vec,
312 const Eigen::VectorXd& scaled_primal_solution,
313 const Eigen::VectorXd& scaled_dual_solution, PointType candidate_type) {
314 const QuadraticProgram& qp = scaled_sharded_qp.Qp();
315 CHECK_EQ(col_scaling_vec.size(), scaled_sharded_qp.PrimalSize());
316 CHECK_EQ(row_scaling_vec.size(), scaled_sharded_qp.DualSize());
317 CHECK_EQ(scaled_primal_solution.size(), scaled_sharded_qp.PrimalSize());
318 CHECK_EQ(scaled_dual_solution.size(), scaled_sharded_qp.DualSize());
319
320 // See https://developers.google.com/optimization/lp/pdlp_math#rescaling for
321 // notes describing the connection between the scaled and unscaled problem.
322
323 ConvergenceInformation result;
324 ResidualNorms primal_residuals = PrimalResidualNorms(
325 scaled_sharded_qp, row_scaling_vec, scaled_primal_solution);
326 result.set_l_inf_primal_residual(primal_residuals.l_inf_residual);
327 result.set_l2_primal_residual(primal_residuals.l_2_residual);
328
329 result.set_l_inf_primal_variable(
330 ScaledLInfNorm(scaled_primal_solution, col_scaling_vec,
331 scaled_sharded_qp.PrimalSharder()));
332 result.set_l2_primal_variable(ScaledNorm(scaled_primal_solution,
333 col_scaling_vec,
334 scaled_sharded_qp.PrimalSharder()));
335 result.set_l_inf_dual_variable(ScaledLInfNorm(
336 scaled_dual_solution, row_scaling_vec, scaled_sharded_qp.DualSharder()));
337 result.set_l2_dual_variable(ScaledNorm(scaled_dual_solution, row_scaling_vec,
338 scaled_sharded_qp.DualSharder()));
339
340 VectorXd scaled_objective_product =
341 ObjectiveProduct(scaled_sharded_qp, scaled_primal_solution);
342 const double quadratic_objective = QuadraticObjective(
343 scaled_sharded_qp, scaled_primal_solution, scaled_objective_product);
344 VectorXd scaled_primal_gradient = PrimalGradientFromObjectiveProduct(
345 scaled_sharded_qp, scaled_dual_solution,
346 std::move(scaled_objective_product));
347 result.set_primal_objective(qp.ApplyObjectiveScalingAndOffset(
348 quadratic_objective + Dot(qp.objective_vector, scaled_primal_solution,
349 scaled_sharded_qp.PrimalSharder())));
350
351 // This is the dual objective from
352 // https://developers.google.com/optimization/lp/pdlp_math minus the last term
353 // (involving r). All scaling terms cancel out.
354 const double dual_objective_piece =
355 -quadratic_objective +
356 DualObjectiveBoundsTerm(scaled_sharded_qp, scaled_dual_solution);
357
358 ResidualNorms dual_residuals =
359 DualResidualNorms(scaled_sharded_qp, col_scaling_vec,
360 scaled_primal_solution, scaled_primal_gradient);
361 result.set_dual_objective(qp.ApplyObjectiveScalingAndOffset(
362 dual_objective_piece + dual_residuals.objective_correction));
363 result.set_corrected_dual_objective(qp.ApplyObjectiveScalingAndOffset(
364 dual_objective_piece + dual_residuals.objective_full_correction));
365 result.set_l_inf_dual_residual(dual_residuals.l_inf_residual);
366 result.set_l2_dual_residual(dual_residuals.l_2_residual);
367
368 result.set_candidate_type(candidate_type);
369 return result;
370}
371
372InfeasibilityInformation ComputeInfeasibilityInformation(
373 const ShardedQuadraticProgram& scaled_sharded_qp,
374 const Eigen::VectorXd& col_scaling_vec,
375 const Eigen::VectorXd& row_scaling_vec,
376 const Eigen::VectorXd& scaled_primal_ray,
377 const Eigen::VectorXd& scaled_dual_ray, PointType candidate_type) {
378 const QuadraticProgram& qp = scaled_sharded_qp.Qp();
379 CHECK_EQ(col_scaling_vec.size(), scaled_sharded_qp.PrimalSize());
380 CHECK_EQ(row_scaling_vec.size(), scaled_sharded_qp.DualSize());
381 CHECK_EQ(scaled_primal_ray.size(), scaled_sharded_qp.PrimalSize());
382 CHECK_EQ(scaled_dual_ray.size(), scaled_sharded_qp.DualSize());
383
384 double l_inf_primal = ScaledLInfNorm(scaled_primal_ray, col_scaling_vec,
385 scaled_sharded_qp.PrimalSharder());
386 double l_inf_dual = ScaledLInfNorm(scaled_dual_ray, row_scaling_vec,
387 scaled_sharded_qp.DualSharder());
388 InfeasibilityInformation result;
389 // Compute primal infeasibility information.
390 VectorXd scaled_primal_gradient = PrimalGradientFromObjectiveProduct(
391 scaled_sharded_qp, scaled_dual_ray,
392 VectorXd::Zero(scaled_sharded_qp.PrimalSize()),
393 /*use_zero_primal_objective=*/true);
394 ResidualNorms dual_residuals =
395 DualResidualNorms(scaled_sharded_qp, col_scaling_vec, scaled_primal_ray,
396 scaled_primal_gradient);
397
398 double dual_ray_objective =
399 DualObjectiveBoundsTerm(scaled_sharded_qp, scaled_dual_ray) +
400 dual_residuals.objective_correction;
401 if (l_inf_dual > 0) {
402 result.set_dual_ray_objective(dual_ray_objective / l_inf_dual);
403 result.set_max_dual_ray_infeasibility(dual_residuals.l_inf_residual /
404 l_inf_dual);
405 } else {
406 result.set_dual_ray_objective(0.0);
407 result.set_max_dual_ray_infeasibility(0.0);
408 }
409
410 // Compute dual infeasibility information.
411 ResidualNorms primal_residuals =
412 PrimalResidualNorms(scaled_sharded_qp, row_scaling_vec, scaled_primal_ray,
413 /*use_homogeneous_constraint_bounds=*/true);
414 // primal_residuals contains the violations of the linear constraints. The
415 // signs of the components are also constrained by the presence or absence
416 // of variable bounds.
417 VectorXd primal_ray_local_sign_max_violation(
418 scaled_sharded_qp.PrimalSharder().NumShards());
419 scaled_sharded_qp.PrimalSharder().ParallelForEachShard(
420 [&](const Sharder::Shard& shard) {
421 const auto lower_bound_shard =
422 shard(scaled_sharded_qp.Qp().variable_lower_bounds);
423 const auto upper_bound_shard =
424 shard(scaled_sharded_qp.Qp().variable_upper_bounds);
425 const auto ray_shard = shard(scaled_primal_ray);
426 const auto scale_shard = shard(col_scaling_vec);
427 double local_max = 0.0;
428 for (int64_t i = 0; i < ray_shard.size(); ++i) {
429 if (std::isfinite(lower_bound_shard[i])) {
430 local_max = std::max(local_max, -ray_shard[i] * scale_shard[i]);
431 }
432 if (std::isfinite(upper_bound_shard[i])) {
433 local_max = std::max(local_max, ray_shard[i] * scale_shard[i]);
434 }
435 }
436 primal_ray_local_sign_max_violation[shard.Index()] = local_max;
437 });
438 const double primal_ray_sign_max_violation =
439 primal_ray_local_sign_max_violation.lpNorm<Eigen::Infinity>();
440
441 if (l_inf_primal > 0.0) {
442 VectorXd scaled_objective_product =
443 ObjectiveProduct(scaled_sharded_qp, scaled_primal_ray);
444 result.set_primal_ray_quadratic_norm(
445 LInfNorm(scaled_objective_product, scaled_sharded_qp.PrimalSharder()) /
446 l_inf_primal);
447 result.set_max_primal_ray_infeasibility(
448 std::max(primal_residuals.l_inf_residual,
449 primal_ray_sign_max_violation) /
450 l_inf_primal);
451 result.set_primal_ray_linear_objective(
452 Dot(scaled_primal_ray, qp.objective_vector,
453 scaled_sharded_qp.PrimalSharder()) /
454 l_inf_primal);
455 } else {
456 result.set_primal_ray_quadratic_norm(0.0);
457 result.set_max_primal_ray_infeasibility(0.0);
458 result.set_primal_ray_linear_objective(0.0);
459 }
460
461 result.set_candidate_type(candidate_type);
462 return result;
463}
464
466 const ShardedQuadraticProgram& sharded_qp, const VectorXd& primal_solution,
467 const VectorXd& dual_solution, PointType candidate_type) {
469 sharded_qp, VectorXd::Ones(sharded_qp.PrimalSize()),
470 VectorXd::Ones(sharded_qp.DualSize()), primal_solution, dual_solution,
471 candidate_type);
472}
473
474VectorXd ReducedCosts(const ShardedQuadraticProgram& sharded_qp,
475 const VectorXd& primal_solution,
476 const VectorXd& dual_solution,
477 bool use_zero_primal_objective) {
478 VectorXd objective_product;
479 if (use_zero_primal_objective) {
480 objective_product = VectorXd::Zero(sharded_qp.PrimalSize());
481 } else {
482 objective_product = ObjectiveProduct(sharded_qp, primal_solution);
483 }
484 VectorXd reduced_costs = PrimalGradientFromObjectiveProduct(
485 sharded_qp, dual_solution, std::move(objective_product),
486 use_zero_primal_objective);
488 [&](const Sharder::Shard& shard) {
489 auto rc_shard = shard(reduced_costs);
490 const auto lower_bound_shard =
491 shard(sharded_qp.Qp().variable_lower_bounds);
492 const auto upper_bound_shard =
493 shard(sharded_qp.Qp().variable_upper_bounds);
494 const auto primal_solution_shard = shard(primal_solution);
495 for (int64_t i = 0; i < rc_shard.size(); ++i) {
496 if (rc_shard[i] != 0.0 &&
497 !HandlePrimalGradientTermAsReducedCost(
498 rc_shard[i], primal_solution_shard[i], lower_bound_shard[i],
499 upper_bound_shard[i])) {
500 rc_shard[i] = 0.0;
501 }
502 }
503 });
504 return reduced_costs;
505}
506
507absl::optional<ConvergenceInformation> GetConvergenceInformation(
508 const IterationStats& stats, PointType candidate_type) {
509 for (const auto& convergence_information : stats.convergence_information()) {
510 if (convergence_information.candidate_type() == candidate_type) {
511 return convergence_information;
512 }
513 }
514 return absl::nullopt;
515}
516
517absl::optional<InfeasibilityInformation> GetInfeasibilityInformation(
518 const IterationStats& stats, PointType candidate_type) {
519 for (const auto& infeasibility_information :
520 stats.infeasibility_information()) {
521 if (infeasibility_information.candidate_type() == candidate_type) {
522 return infeasibility_information;
523 }
524 }
525 return absl::nullopt;
526}
527
528absl::optional<PointMetadata> GetPointMetadata(const IterationStats& stats,
529 const PointType point_type) {
530 for (const auto& metadata : stats.point_metadata()) {
531 if (metadata.point_type() == point_type) {
532 return metadata;
533 }
534 }
535 return absl::nullopt;
536}
537
539 const Eigen::VectorXd& primal_solution,
540 const Eigen::VectorXd& dual_solution,
541 const std::vector<int>& random_projection_seeds,
542 PointMetadata& metadata) {
543 for (const int random_projection_seed : random_projection_seeds) {
544 std::mt19937 seed_generator(random_projection_seed);
545 metadata.mutable_random_primal_projections()->Add(RandomProjection(
546 primal_solution, sharded_qp.PrimalSharder(), seed_generator));
547 metadata.mutable_random_dual_projections()->Add(RandomProjection(
548 dual_solution, sharded_qp.DualSharder(), seed_generator));
549 }
550}
551
552} // namespace operations_research::pdlp
int64_t max
Definition: alldiff_cst.cc:140
#define CHECK_EQ(val1, val2)
Definition: base/logging.h:703
static T Square(const T x)
Definition: mathutil.h:101
void ParallelForEachShard(const std::function< void(const Shard &)> &func) const
Definition: sharder.cc:97
double upper_bound
double lower_bound
double l_2_residual
double objective_full_correction
double objective_correction
double l_inf_residual
double ScaledNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition: sharder.cc:253
double Dot(const VectorXd &v1, const VectorXd &v2, const Sharder &sharder)
Definition: sharder.cc:197
ConvergenceInformation ComputeScaledConvergenceInformation(const ShardedQuadraticProgram &sharded_qp, const VectorXd &primal_solution, const VectorXd &dual_solution, PointType candidate_type)
double LInfNorm(const VectorXd &vector, const Sharder &sharder)
Definition: sharder.cc:202
VectorXd TransposedMatrixVectorProduct(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &vector, const Sharder &sharder)
Definition: sharder.cc:151
double ScaledLInfNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition: sharder.cc:236
VectorXd ReducedCosts(const ShardedQuadraticProgram &sharded_qp, const VectorXd &primal_solution, const VectorXd &dual_solution, bool use_zero_primal_objective)
absl::optional< PointMetadata > GetPointMetadata(const IterationStats &stats, const PointType point_type)
InfeasibilityInformation ComputeInfeasibilityInformation(const ShardedQuadraticProgram &scaled_sharded_qp, const Eigen::VectorXd &col_scaling_vec, const Eigen::VectorXd &row_scaling_vec, const Eigen::VectorXd &scaled_primal_ray, const Eigen::VectorXd &scaled_dual_ray, PointType candidate_type)
absl::optional< ConvergenceInformation > GetConvergenceInformation(const IterationStats &stats, PointType candidate_type)
bool IsLinearProgram(const QuadraticProgram &qp)
void SetRandomProjections(const ShardedQuadraticProgram &sharded_qp, const Eigen::VectorXd &primal_solution, const Eigen::VectorXd &dual_solution, const std::vector< int > &random_projection_seeds, PointMetadata &metadata)
absl::optional< InfeasibilityInformation > GetInfeasibilityInformation(const IterationStats &stats, PointType candidate_type)
ConvergenceInformation ComputeConvergenceInformation(const ShardedQuadraticProgram &scaled_sharded_qp, const Eigen::VectorXd &col_scaling_vec, const Eigen::VectorXd &row_scaling_vec, const Eigen::VectorXd &scaled_primal_solution, const Eigen::VectorXd &scaled_dual_solution, PointType candidate_type)
int64_t Zero()
NOLINT.
double ApplyObjectiveScalingAndOffset(double objective) const