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"
43using ::Eigen::VectorXd;
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());
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);
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]))
91 : upper_bound_shard[i];
92 const double lower_bound = (use_homogeneous_constraint_bounds &&
93 std::isfinite(lower_bound_shard[i]))
95 : lower_bound_shard[i];
96 double scaled_residual = 0.0;
98 scaled_residual = primal_product_shard[i] -
upper_bound;
100 }
else if (primal_product_shard[i] <
lower_bound) {
101 scaled_residual =
lower_bound - primal_product_shard[i];
103 const double residual = scaled_residual / row_scaling_shard[i];
105 sumsq_residual += residual * residual;
108 local_sumsq_residual[shard.Index()] = sumsq_residual;
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>(),
120bool HandlePrimalGradientTermAsReducedCost(
double primal_gradient,
124 if (primal_gradient == 0.0)
return true;
125 return std::abs(primal_value -
127 std::abs(primal_value);
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;
162 double sumsq_residual = 0.0;
163 for (int64_t i = 0; i < primal_gradient_shard.size(); ++i) {
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];
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;
184 local_dual_correction[shard.Index()] = dual_correction;
185 local_dual_full_correction[shard.Index()] = dual_full_correction;
187 local_sumsq_residual[shard.Index()] = sumsq_residual;
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>(),
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());
205 sharded_qp.PrimalSharder().ParallelForEachShard(
206 [&](
const Sharder::Shard& shard) {
208 shard(*sharded_qp.Qp().objective_matrix) * shard(primal_solution);
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());
221 Dot(objective_product, primal_solution, sharded_qp.PrimalSharder());
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());
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;
242 shard(objective_product) +=
243 shard(qp.objective_vector) -
244 shard(qp.constraint_matrix).transpose() * dual_solution;
247 return objective_product;
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) {
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);
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];
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)());
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];
301 dot_product[shard.Index()] = shard_dot_product;
302 gaussian_norm_squared[shard.Index()] = shard_norm_squared;
304 return dot_product.sum() / std::sqrt(gaussian_norm_squared.sum());
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) {
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);
329 result.set_l_inf_primal_variable(
332 result.set_l2_primal_variable(
ScaledNorm(scaled_primal_solution,
336 scaled_dual_solution, row_scaling_vec, scaled_sharded_qp.
DualSharder()));
337 result.set_l2_dual_variable(
ScaledNorm(scaled_dual_solution, row_scaling_vec,
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));
354 const double dual_objective_piece =
355 -quadratic_objective +
356 DualObjectiveBoundsTerm(scaled_sharded_qp, scaled_dual_solution);
358 ResidualNorms dual_residuals =
359 DualResidualNorms(scaled_sharded_qp, col_scaling_vec,
360 scaled_primal_solution, scaled_primal_gradient);
362 dual_objective_piece + dual_residuals.objective_correction));
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);
368 result.set_candidate_type(candidate_type);
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) {
384 double l_inf_primal =
ScaledLInfNorm(scaled_primal_ray, col_scaling_vec,
386 double l_inf_dual =
ScaledLInfNorm(scaled_dual_ray, row_scaling_vec,
388 InfeasibilityInformation result;
390 VectorXd scaled_primal_gradient = PrimalGradientFromObjectiveProduct(
391 scaled_sharded_qp, scaled_dual_ray,
394 ResidualNorms dual_residuals =
395 DualResidualNorms(scaled_sharded_qp, col_scaling_vec, scaled_primal_ray,
396 scaled_primal_gradient);
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 /
406 result.set_dual_ray_objective(0.0);
407 result.set_max_dual_ray_infeasibility(0.0);
411 ResidualNorms primal_residuals =
412 PrimalResidualNorms(scaled_sharded_qp, row_scaling_vec, scaled_primal_ray,
417 VectorXd primal_ray_local_sign_max_violation(
421 const auto lower_bound_shard =
423 const auto upper_bound_shard =
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]);
432 if (std::isfinite(upper_bound_shard[i])) {
433 local_max =
std::max(local_max, ray_shard[i] * scale_shard[i]);
436 primal_ray_local_sign_max_violation[shard.
Index()] = local_max;
438 const double primal_ray_sign_max_violation =
439 primal_ray_local_sign_max_violation.lpNorm<Eigen::Infinity>();
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(
447 result.set_max_primal_ray_infeasibility(
448 std::max(primal_residuals.l_inf_residual,
449 primal_ray_sign_max_violation) /
451 result.set_primal_ray_linear_objective(
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);
461 result.set_candidate_type(candidate_type);
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,
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) {
482 objective_product = ObjectiveProduct(sharded_qp, primal_solution);
484 VectorXd reduced_costs = PrimalGradientFromObjectiveProduct(
485 sharded_qp, dual_solution, std::move(objective_product),
486 use_zero_primal_objective);
489 auto rc_shard = shard(reduced_costs);
490 const auto lower_bound_shard =
492 const auto upper_bound_shard =
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])) {
504 return reduced_costs;
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;
514 return absl::nullopt;
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;
525 return absl::nullopt;
529 const PointType point_type) {
530 for (
const auto& metadata : stats.point_metadata()) {
531 if (metadata.point_type() == point_type) {
535 return absl::nullopt;
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));
#define CHECK_EQ(val1, val2)
static T Square(const T x)
const Sharder & DualSharder() const
const Sharder & PrimalSharder() const
int64_t PrimalSize() const
const QuadraticProgram & Qp() const
void ParallelForEachShard(const std::function< void(const Shard &)> &func) const
double objective_full_correction
double objective_correction
double ScaledNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
double Dot(const VectorXd &v1, const VectorXd &v2, const Sharder &sharder)
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)
VectorXd TransposedMatrixVectorProduct(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &vector, const Sharder &sharder)
double ScaledLInfNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
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)
Eigen::VectorXd variable_upper_bounds
Eigen::VectorXd variable_lower_bounds
double ApplyObjectiveScalingAndOffset(double objective) const
Eigen::VectorXd objective_vector