23#include "Eigen/SparseCore"
24#include "absl/types/optional.h"
25#include "gmock/gmock.h"
26#include "gtest/gtest.h"
30#include "ortools/pdlp/solve_log.pb.h"
36using ::Eigen::VectorXd;
37using ::testing::ElementsAre;
38using ::testing::ElementsAreArray;
39using ::testing::IsNan;
41TEST(ShardedWeightedAverageTest, SimpleAverage) {
44 Eigen::VectorXd vec1(2), vec2(2);
48 ShardedWeightedAverage
average(&sharder);
52 ASSERT_TRUE(
average.HasNonzeroWeight());
53 EXPECT_EQ(
average.NumTerms(), 2);
55 EXPECT_THAT(
average.ComputeAverage(), ElementsAre(2.0, 5.0));
58 EXPECT_FALSE(
average.HasNonzeroWeight());
59 EXPECT_EQ(
average.NumTerms(), 0);
62TEST(ShardedWeightedAverageTest, MoveConstruction) {
65 Eigen::VectorXd vec(2);
68 ShardedWeightedAverage
average(&sharder);
71 ShardedWeightedAverage average2(std::move(
average));
72 EXPECT_THAT(average2.ComputeAverage(), ElementsAre(4.0, 1.0));
75TEST(ShardedWeightedAverageTest, MoveAssignment) {
76 Sharder sharder1(2, 2,
78 Sharder sharder2(3, 2,
80 Eigen::VectorXd vec1(2), vec2(2);
84 ShardedWeightedAverage average1(&sharder1);
85 average1.Add(vec1, 2.0);
87 ShardedWeightedAverage average2(&sharder2);
89 average2 = std::move(average1);
90 average2.Add(vec2, 2.0);
91 EXPECT_THAT(average2.ComputeAverage(), ElementsAre(2.0, 2.0));
94TEST(ShardedWeightedAverageTest, ZeroAverage) {
98 ShardedWeightedAverage
average(&sharder);
99 ASSERT_FALSE(
average.HasNonzeroWeight());
101 EXPECT_THAT(
average.ComputeAverage(), ElementsAre(0.0));
106TEST(ShardedWeightedAverageTest, AveragesEqualWithoutRoundoff) {
107 Sharder sharder(4, 1,
109 ShardedWeightedAverage
average(&sharder);
110 EXPECT_THAT(
average.ComputeAverage(), ElementsAre(0, 0, 0, 0));
112 data << 1.0, 1.0 / 3, 3.0 / 7, 3.14159;
114 EXPECT_THAT(
average.ComputeAverage(), ElementsAreArray(data));
116 EXPECT_THAT(
average.ComputeAverage(), ElementsAreArray(data));
118 EXPECT_THAT(
average.ComputeAverage(), ElementsAreArray(data));
126 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
129 EXPECT_EQ(stats.num_variables(), 4);
130 EXPECT_EQ(stats.num_constraints(), 4);
131 EXPECT_DOUBLE_EQ(stats.constraint_matrix_col_min_l_inf_norm(), 1.0);
132 EXPECT_DOUBLE_EQ(stats.constraint_matrix_row_min_l_inf_norm(), 1.0);
133 EXPECT_EQ(stats.constraint_matrix_num_nonzeros(), 9);
134 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_max(), 4.0);
135 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_min(), 1.0);
136 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_avg(), 14.5 / 9.0);
137 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_max(), 5.5);
138 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_min(), 1.0);
139 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_avg(), 2.375);
140 EXPECT_DOUBLE_EQ(stats.objective_vector_l2_norm(), std::sqrt(36.25));
141 EXPECT_EQ(stats.objective_matrix_num_nonzeros(), 0);
142 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_max(), 0.0);
143 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_min(), 0.0);
144 EXPECT_THAT(stats.objective_matrix_abs_avg(), IsNan());
145 EXPECT_EQ(stats.variable_bound_gaps_num_finite(), 1);
146 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_max(), 1.0);
147 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_min(), 1.0);
148 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_avg(), 1.0);
149 EXPECT_DOUBLE_EQ(stats.combined_bounds_max(), 12.0);
150 EXPECT_DOUBLE_EQ(stats.combined_bounds_min(), 1.0);
151 EXPECT_DOUBLE_EQ(stats.combined_bounds_avg(), 6.0);
152 EXPECT_DOUBLE_EQ(stats.combined_bounds_l2_norm(), std::sqrt(210.0));
156 ShardedQuadraticProgram lp(
TinyLp(), 2, 2);
159 EXPECT_EQ(stats.num_variables(), 4);
160 EXPECT_EQ(stats.num_constraints(), 3);
161 EXPECT_DOUBLE_EQ(stats.constraint_matrix_col_min_l_inf_norm(), 1.0);
162 EXPECT_DOUBLE_EQ(stats.constraint_matrix_row_min_l_inf_norm(), 1.0);
163 EXPECT_EQ(stats.constraint_matrix_num_nonzeros(), 8);
164 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_max(), 2.0);
165 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_min(), 1.0);
166 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_avg(), 1.25);
167 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_max(), 5.0);
168 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_min(), 1.0);
169 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_avg(), 2.25);
170 EXPECT_DOUBLE_EQ(stats.objective_vector_l2_norm(), std::sqrt(31.0));
171 EXPECT_EQ(stats.objective_matrix_num_nonzeros(), 0);
172 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_max(), 0.0);
173 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_min(), 0.0);
174 EXPECT_THAT(stats.objective_matrix_abs_avg(), IsNan());
175 EXPECT_EQ(stats.variable_bound_gaps_num_finite(), 4);
176 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_max(), 6.0);
177 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_min(), 2.0);
178 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_avg(), 3.75);
179 EXPECT_DOUBLE_EQ(stats.combined_bounds_max(), 12.0);
180 EXPECT_DOUBLE_EQ(stats.combined_bounds_min(), 1.0);
181 EXPECT_DOUBLE_EQ(stats.combined_bounds_avg(), 20.0 / 3.0);
182 EXPECT_DOUBLE_EQ(stats.combined_bounds_l2_norm(), std::sqrt(194.0));
189 EXPECT_EQ(stats.num_variables(), 2);
190 EXPECT_EQ(stats.num_constraints(), 1);
191 EXPECT_DOUBLE_EQ(stats.constraint_matrix_col_min_l_inf_norm(), 1.0);
192 EXPECT_DOUBLE_EQ(stats.constraint_matrix_row_min_l_inf_norm(), 1.0);
193 EXPECT_EQ(stats.constraint_matrix_num_nonzeros(), 2);
194 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_max(), 1.0);
195 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_min(), 1.0);
196 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_avg(), 1.0);
197 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_max(), 1.0);
198 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_min(), 1.0);
199 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_avg(), 1.0);
200 EXPECT_DOUBLE_EQ(stats.objective_vector_l2_norm(), std::sqrt(2.0));
201 EXPECT_EQ(stats.objective_matrix_num_nonzeros(), 2);
202 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_max(), 4.0);
203 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_min(), 1.0);
204 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_avg(), 2.5);
205 EXPECT_EQ(stats.variable_bound_gaps_num_finite(), 2);
206 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_max(), 6.0);
207 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_min(), 1.0);
208 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_avg(), 3.5);
209 EXPECT_DOUBLE_EQ(stats.combined_bounds_max(), 1.0);
210 EXPECT_DOUBLE_EQ(stats.combined_bounds_min(), 1.0);
211 EXPECT_DOUBLE_EQ(stats.combined_bounds_avg(), 1.0);
212 EXPECT_DOUBLE_EQ(stats.combined_bounds_l2_norm(), 1.0);
215TEST(ProblemStatsTest, ModifiedTestDiagonalQp1) {
218 orig_qp.objective_matrix->diagonal() << 2.0, 0.0;
219 ShardedQuadraticProgram qp(orig_qp, 2, 2);
222 EXPECT_EQ(stats.objective_matrix_num_nonzeros(), 1);
223 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_max(), 2.0);
224 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_min(), 0.0);
225 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_avg(), 1.0);
230TEST(ProblemStatsTest, TestLpWithInfiniteConstraintBoundThreshold) {
231 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
232 const QuadraticProgramStats stats =
235 EXPECT_EQ(stats.num_variables(), 4);
236 EXPECT_EQ(stats.num_constraints(), 4);
237 EXPECT_DOUBLE_EQ(stats.constraint_matrix_col_min_l_inf_norm(), 1.0);
238 EXPECT_DOUBLE_EQ(stats.constraint_matrix_row_min_l_inf_norm(), 1.0);
239 EXPECT_EQ(stats.constraint_matrix_num_nonzeros(), 9);
240 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_max(), 4.0);
241 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_min(), 1.0);
242 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_avg(), 14.5 / 9.0);
243 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_max(), 5.5);
244 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_min(), 1.0);
245 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_avg(), 2.375);
246 EXPECT_DOUBLE_EQ(stats.objective_vector_l2_norm(), std::sqrt(36.25));
247 EXPECT_EQ(stats.objective_matrix_num_nonzeros(), 0);
248 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_max(), 0.0);
249 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_min(), 0.0);
250 EXPECT_THAT(stats.objective_matrix_abs_avg(), IsNan());
251 EXPECT_EQ(stats.variable_bound_gaps_num_finite(), 1);
252 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_max(), 1.0);
253 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_min(), 1.0);
254 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_avg(), 1.0);
255 EXPECT_DOUBLE_EQ(stats.combined_bounds_max(), 7.0);
256 EXPECT_DOUBLE_EQ(stats.combined_bounds_min(), 0.0);
257 EXPECT_DOUBLE_EQ(stats.combined_bounds_avg(), 3.0);
258 EXPECT_DOUBLE_EQ(stats.combined_bounds_l2_norm(), std::sqrt(66.0));
261TEST(ProblemStatsTest, NoFiniteGaps) {
265 EXPECT_EQ(stats.variable_bound_gaps_num_finite(), 0);
266 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_max(), 0.0);
267 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_min(), 0.0);
268 EXPECT_THAT(stats.variable_bound_gaps_avg(), IsNan());
276 EXPECT_EQ(stats.constraint_matrix_num_nonzeros(), 0);
277 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_max(), 0.0);
278 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_min(), 0.0);
279 EXPECT_THAT(stats.constraint_matrix_abs_avg(), IsNan());
280 EXPECT_DOUBLE_EQ(stats.constraint_matrix_col_min_l_inf_norm(), 0.0);
281 EXPECT_DOUBLE_EQ(stats.constraint_matrix_row_min_l_inf_norm(), 0.0);
282 EXPECT_DOUBLE_EQ(stats.combined_bounds_max(), 0.0);
283 EXPECT_DOUBLE_EQ(stats.combined_bounds_min(), 0.0);
284 EXPECT_THAT(stats.combined_bounds_avg(), IsNan());
287TEST(ProblemStatsTest, EmptyLp) {
288 ShardedQuadraticProgram lp(QuadraticProgram(0, 0), 2, 2);
292 EXPECT_EQ(stats.num_variables(), 0);
293 EXPECT_EQ(stats.num_constraints(), 0);
294 EXPECT_DOUBLE_EQ(stats.constraint_matrix_col_min_l_inf_norm(), 0.0);
295 EXPECT_DOUBLE_EQ(stats.constraint_matrix_row_min_l_inf_norm(), 0.0);
296 EXPECT_EQ(stats.constraint_matrix_num_nonzeros(), 0);
297 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_max(), 0.0);
298 EXPECT_DOUBLE_EQ(stats.constraint_matrix_abs_min(), 0.0);
299 EXPECT_THAT(stats.constraint_matrix_abs_avg(), IsNan());
300 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_max(), 0.0);
301 EXPECT_DOUBLE_EQ(stats.objective_vector_abs_min(), 0.0);
302 EXPECT_THAT(stats.objective_vector_abs_avg(), IsNan());
303 EXPECT_DOUBLE_EQ(stats.objective_vector_l2_norm(), 0.0);
304 EXPECT_EQ(stats.objective_matrix_num_nonzeros(), 0);
305 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_max(), 0.0);
306 EXPECT_DOUBLE_EQ(stats.objective_matrix_abs_min(), 0.0);
307 EXPECT_THAT(stats.objective_matrix_abs_avg(), IsNan());
308 EXPECT_EQ(stats.variable_bound_gaps_num_finite(), 0);
309 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_max(), 0.0);
310 EXPECT_DOUBLE_EQ(stats.variable_bound_gaps_min(), 0.0);
311 EXPECT_THAT(stats.variable_bound_gaps_avg(), IsNan());
312 EXPECT_DOUBLE_EQ(stats.combined_bounds_max(), 0.0);
313 EXPECT_DOUBLE_EQ(stats.combined_bounds_min(), 0.0);
314 EXPECT_THAT(stats.combined_bounds_avg(), IsNan());
315 EXPECT_DOUBLE_EQ(stats.combined_bounds_l2_norm(), 0.0);
323 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
324 VectorXd row_scaling_vec(4), col_scaling_vec(4);
325 row_scaling_vec << 1, 2, 1, 3;
326 col_scaling_vec << 0, 1, 2, -1;
328 EXPECT_THAT(row_scaling_vec, ElementsAre(1 / std::sqrt(2), 1.0, 1.0, 1.0));
329 EXPECT_THAT(col_scaling_vec,
330 ElementsAre(0.0, 1.0, 2.0 / 3.0, -1.0 / std::sqrt(3.0)));
337TEST(L2RuizRescaling, OneIteration) {
338 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
339 VectorXd row_scaling_vec(4), col_scaling_vec(4);
340 row_scaling_vec << 1, 2, 1, 3;
341 col_scaling_vec << 0, 1, 2, -1;
343 EXPECT_THAT(row_scaling_vec, ElementsAre(1.0 / std::pow(3.0, 0.5), 1.0, 1.0,
344 3.0 / std::pow(90.0, 0.25)));
345 EXPECT_THAT(col_scaling_vec, ElementsAre(0.0, 1.0, 2.0 / std::pow(101, 0.25),
346 -1.0 / std::pow(13.0, 0.25)));
351TEST(L2RuizRescaling, OneIterationNonSquare) {
352 QuadraticProgram test_lp(2, 1);
353 std::vector<Eigen::Triplet<double, int64_t>> triplets = {{0, 0, 2.0},
355 test_lp.constraint_matrix.setFromTriplets(triplets.begin(), triplets.end());
356 ShardedQuadraticProgram lp(std::move(test_lp), 2,
358 VectorXd row_scaling_vec = VectorXd::Ones(1);
359 VectorXd col_scaling_vec = VectorXd::Ones(2);
361 EXPECT_THAT(row_scaling_vec, ElementsAre(1.0 / std::pow(13.0, 0.25)));
362 EXPECT_THAT(col_scaling_vec,
363 ElementsAre(1.0 / std::sqrt(2.0), 1.0 / std::sqrt(3.0)));
369 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
370 VectorXd row_scaling_vec(4), col_scaling_vec(4);
371 VectorXd col_norm(4), row_norm(4);
372 row_scaling_vec << 1, 1, 1, 1;
373 col_scaling_vec << 1, 1, 1, 1;
377 col_scaling_vec, lp.ConstraintMatrixSharder());
380 lp.TransposedConstraintMatrixSharder());
381 EXPECT_THAT(row_norm, EigenArrayNear<double>({1.0, 1.0, 1.0, 1.0}, 1.0e-4));
382 EXPECT_THAT(col_norm, EigenArrayNear<double>({1.0, 1.0, 1.0, 1.0}, 1.0e-4));
394 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
396 RescalingOptions{.l_inf_ruiz_iterations = 1, .l2_norm_rescaling =
true},
398 EXPECT_THAT(scaling.row_scaling_vec,
399 EigenArrayNear<double>(
400 {1.0 / sqrt(2.0 * 1.5275), 1.0 / sqrt(1.0 * 0.9574),
401 1.0 / sqrt(4.0 * 1.0), 1.0 / sqrt(1.5 * 1.1547)},
403 EXPECT_THAT(scaling.col_scaling_vec,
404 EigenArrayNear<double>(
405 {1.0 / sqrt(4.0 * 1.3229), 1.0 / sqrt(1.0 * 0.7071),
406 1.0 / sqrt(1.5 * 1.4142), 1.0 / sqrt(2.0 * 1.1547)},
410TEST(ComputePrimalGradientTest, CorrectForLp) {
413 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
415 VectorXd primal_solution(4), dual_solution(4);
416 primal_solution << 0.0, 0.0, 0.0, 3.0;
417 dual_solution << -1.0, 0.0, 1.0, 1.0;
420 lp, primal_solution, lp.TransposedConstraintMatrix() * dual_solution);
424 EXPECT_THAT(primal_part.gradient,
425 ElementsAre(5.5 - 2.0, -2.0 + 1.0, -1.0 - 0.5, 1.0 + 3.0));
427 EXPECT_DOUBLE_EQ(primal_part.value, 3.0 + 9.0);
430TEST(ComputeDualGradientTest, CorrectForLp) {
433 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
435 VectorXd primal_solution(4), dual_solution(4);
436 primal_solution << 0.0, 0.0, 0.0, 3.0;
437 dual_solution << -1.0, 0.0, 1.0, 1.0;
440 lp, dual_solution, lp.Qp().constraint_matrix * primal_solution);
444 EXPECT_THAT(dual_part.gradient,
445 ElementsAre(12.0 - 6.0, 7.0, -4.0, -1.0 + 3.0));
447 EXPECT_DOUBLE_EQ(dual_part.value, 12.0 * -1.0 + -4.0 * 1.0 + -1.0 * 1.0);
450TEST(ComputeDualGradientTest, CorrectOnTwoSidedConstraints) {
451 QuadraticProgram qp =
TestLp();
455 qp.constraint_lower_bounds[0] = 4;
456 qp.constraint_lower_bounds[1] = 5;
457 qp.constraint_upper_bounds[2] = -1;
458 ShardedQuadraticProgram sharded_qp(std::move(qp), 2,
461 VectorXd primal_solution(4), dual_solution(4);
462 primal_solution << 0.0, 0.0, 0.0, 3.0;
463 dual_solution << 0.0, 0.0, 0.0, -1.0;
465 const LagrangianPart dual_part =
467 sharded_qp.Qp().constraint_matrix * primal_solution);
471 EXPECT_THAT(dual_part.gradient,
472 ElementsAre(0.0, 5.0 - 0.0, -1.0 - 0.0, 1.0 + 3.0));
474 EXPECT_DOUBLE_EQ(dual_part.value, 1.0 * -1.0);
477TEST(HasValidBoundsTest, SmallInvalidLp) {
482 EXPECT_FALSE(is_valid);
485TEST(HasValidBoundsTest, SmallValidLp) {
490 EXPECT_TRUE(is_valid);
493TEST(ComputePrimalGradientTest, CorrectForQp) {
497 VectorXd primal_solution(2), dual_solution(1);
498 primal_solution << 1.0, 2.0;
499 dual_solution << -2.0;
502 qp, primal_solution, qp.TransposedConstraintMatrix() * dual_solution);
507 EXPECT_THAT(primal_part.gradient,
508 ElementsAre(-1.0 + 2.0 + 4.0, -1.0 + 2.0 + 2.0));
510 EXPECT_DOUBLE_EQ(primal_part.value, 4.0 - 3.0 + 2.0 * 3.0);
513TEST(ComputeDualGradientTest, CorrectForQp) {
517 VectorXd primal_solution(2), dual_solution(1);
518 primal_solution << 1.0, 2.0;
519 dual_solution << -2.0;
522 qp, dual_solution, qp.Qp().constraint_matrix * primal_solution);
527 EXPECT_THAT(dual_part.gradient, ElementsAre(1.0 - (1.0 + 2.0)));
529 EXPECT_DOUBLE_EQ(dual_part.value, -2.0);
532TEST(EstimateSingularValuesTest, CorrectForTestLp) {
533 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
536 std::mt19937 random(1);
538 lp, absl::nullopt, absl::nullopt,
541 EXPECT_NEAR(result.singular_value, 4.76945, 0.01);
542 EXPECT_LT(result.num_iterations, 300);
545TEST(EstimateSingularValuesTest, CorrectForTestLpWithActivePrimalSubspace) {
546 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
548 VectorXd primal_solution(4);
550 primal_solution << 0.0, -2.0, 0.0, 3.0;
553 std::mt19937 random(1);
555 lp, primal_solution, absl::nullopt, 0.01,
557 EXPECT_NEAR(result.singular_value, 4.73818, 0.01);
558 EXPECT_LT(result.num_iterations, 300);
561TEST(EstimateSingularValuesTest, CorrectForTestLpWithActiveDualSubspace) {
562 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
564 VectorXd dual_solution(4);
567 dual_solution << 1.0, 0.0, 1.0, 3.0;
570 std::mt19937 random(1);
572 lp, absl::nullopt, dual_solution, 0.01,
574 EXPECT_NEAR(result.singular_value, 4.64203, 0.01);
575 EXPECT_LT(result.num_iterations, 300);
578TEST(EstimateSingularValuesTest, CorrectForTestLpWithBothActiveSubspaces) {
579 ShardedQuadraticProgram lp(
TestLp(), 2, 2);
581 VectorXd primal_solution(4), dual_solution(4);
583 primal_solution << 0.0, -2.0, 0.0, 3.0;
586 dual_solution << 1.0, 0.0, 1.0, 3.0;
589 std::mt19937 random(1);
591 lp, primal_solution, dual_solution, 0.01,
593 EXPECT_NEAR(result.singular_value, 4.60829, 0.01);
594 EXPECT_LT(result.num_iterations, 300);
597TEST(ProjectToPrimalVariableBoundsTest,
TestLp) {
598 ShardedQuadraticProgram qp(
TestLp(), 2,
601 primal << -3, -3, 5, 5;
603 EXPECT_THAT(primal, ElementsAre(-3, -2, 5, 3.5));
606TEST(ProjectToDualVariableBoundsTest,
TestLp) {
607 ShardedQuadraticProgram qp(
TestLp(), 2,
610 dual << 1, 1, -1, -1;
612 EXPECT_THAT(dual, ElementsAre(1, 0, 0, -1));
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)
VectorXd ScaledColLInfNorm(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &row_scaling_vec, const VectorXd &col_scaling_vec, const Sharder &sharder)
bool HasValidBounds(const QuadraticProgram &qp)
QuadraticProgram TinyLp()
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)
QuadraticProgram LpWithoutConstraints()
QuadraticProgram SmallInvalidProblemLp()
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)
QuadraticProgram TestLp()
QuadraticProgram SmallPrimalInfeasibleLp()
QuadraticProgram TestDiagonalQp1()
TEST(LinearAssignmentTest, NullMatrix)