25#include "absl/log/check.h"
35using ::Eigen::VectorXd;
75class VectorTrustRegionProblem {
77 VectorTrustRegionProblem(
const VectorXd* objective,
80 const VectorXd* center_point,
81 const VectorXd* norm_weight)
82 : objective_(*objective),
85 center_point_(*center_point),
86 norm_weight_(*norm_weight) {}
87 double Objective(int64_t
index)
const {
return objective_(
index); }
90 double CenterPoint(int64_t
index)
const {
return center_point_(
index); }
91 double NormWeight(int64_t
index)
const {
return norm_weight_(
index); }
94 const VectorXd& objective_;
95 const VectorXd& lower_bound_;
96 const VectorXd& upper_bound_;
97 const VectorXd& center_point_;
98 const VectorXd& norm_weight_;
115class JointTrustRegionProblem {
117 JointTrustRegionProblem(
const QuadraticProgram* qp,
119 const VectorXd* dual_solution,
120 const VectorXd* primal_gradient,
121 const VectorXd* dual_gradient,
122 const double primal_weight)
124 primal_size_(qp_.variable_lower_bounds.
size()),
126 dual_solution_(*dual_solution),
127 primal_gradient_(*primal_gradient),
128 dual_gradient_(*dual_gradient),
129 primal_weight_(primal_weight) {}
130 double Objective(int64_t
index)
const {
131 return index < primal_size_ ? primal_gradient_[
index]
132 : -dual_gradient_[
index - primal_size_];
135 return index < primal_size_ ? qp_.variable_lower_bounds[
index]
136 : std::isfinite(qp_.constraint_upper_bounds[
index - primal_size_])
137 ? -std::numeric_limits<double>::infinity()
141 return index < primal_size_ ? qp_.variable_upper_bounds[
index]
142 : std::isfinite(qp_.constraint_lower_bounds[
index - primal_size_])
143 ? std::numeric_limits<double>::infinity()
146 double CenterPoint(int64_t
index)
const {
147 return index < primal_size_ ? primal_solution_[
index]
148 : dual_solution_[
index - primal_size_];
150 double NormWeight(int64_t
index)
const {
151 return index < primal_size_ ? 0.5 * primal_weight_ : 0.5 / primal_weight_;
155 const QuadraticProgram& qp_;
156 const int64_t primal_size_;
157 const VectorXd& primal_solution_;
158 const VectorXd& dual_solution_;
159 const VectorXd& primal_gradient_;
160 const VectorXd& dual_gradient_;
161 const double primal_weight_;
164struct TrustRegionResultStepSize {
177template <
typename TrustRegionProblem>
178double MedianOfShardMedians(
179 const TrustRegionProblem& problem,
180 const std::vector<std::vector<int64_t>>& indexed_components_by_shard,
181 const Sharder& sharder) {
182 std::vector<std::optional<double>> shard_medians(sharder.NumShards(),
184 sharder.ParallelForEachShard([&](
const Sharder::Shard& shard) {
185 const auto& indexed_shard_components =
186 indexed_components_by_shard[shard.Index()];
187 if (!indexed_shard_components.empty()) {
189 indexed_shard_components, [&](
const int64_t
index) {
194 std::vector<double> non_empty_medians;
195 for (
const auto& median : shard_medians) {
196 if (median.has_value()) {
197 non_empty_medians.push_back(*median);
200 CHECK(!non_empty_medians.empty());
202 [](
const double x) {
return x; });
210template <
typename TrustRegionProblem>
211InitialState ComputeInitialState(
const TrustRegionProblem& problem,
212 const Sharder& sharder) {
214 result.undecided_components_by_shard.resize(sharder.NumShards());
215 result.radius_coefficient_of_decided_components =
216 sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
217 const int64_t shard_start = sharder.ShardStart(shard.Index());
218 const int64_t shard_size = sharder.ShardSize(shard.Index());
220 problem, shard_start, shard_start + shard_size,
221 result.undecided_components_by_shard[shard.Index()]);
226template <
typename TrustRegionProblem>
228 const TrustRegionProblem& problem,
const double step_size,
229 const Sharder& sharder,
231 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
237template <
typename TrustRegionProblem>
239 const TrustRegionProblem& problem,
const double step_size_threshold,
240 const Sharder& sharder,
242 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
244 problem, step_size_threshold,
249template <
typename TrustRegionProblem>
251 const TrustRegionProblem& problem,
const double step_size_threshold,
252 const Sharder& sharder,
254 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
256 problem, step_size_threshold,
261int64_t NumUndecidedComponents(
263 int64_t num_undecided_components = 0;
265 num_undecided_components += undecided_components.size();
267 return num_undecided_components;
270int64_t MaxUndecidedComponentsInAnyShard(
274 max = std::max<int64_t>(
max, undecided_components.size());
279template <
typename TrustRegionProblem>
280VectorXd ComputeSolution(
const TrustRegionProblem& problem,
281 const double step_size,
const Sharder& sharder) {
282 VectorXd
solution(sharder.NumElements());
283 sharder.ParallelForEachShard([&](
const Sharder::Shard& shard) {
284 const int64_t shard_start = sharder.ShardStart(shard.Index());
285 const int64_t shard_size = sharder.ShardSize(shard.Index());
286 for (int64_t
index = shard_start;
index < shard_start + shard_size;
294template <
typename TrustRegionProblem>
296 const double step_size,
const Sharder& sharder) {
297 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
298 const int64_t shard_start = sharder.ShardStart(shard.Index());
299 const int64_t shard_size = sharder.ShardSize(shard.Index());
300 double shard_value = 0.0;
301 for (int64_t
index = shard_start;
index < shard_start + shard_size;
303 shard_value += problem.Objective(
index) *
305 problem.CenterPoint(
index));
332template <
typename TrustRegionProblem>
333TrustRegionResultStepSize SolveTrustRegionStepSize(
334 const TrustRegionProblem& problem,
const double target_radius,
335 const Sharder& sharder) {
336 CHECK_GE(target_radius, 0.0);
338 const bool norm_weights_are_positive =
339 sharder.ParallelTrueForAllShards([&](
const Sharder::Shard& shard) {
340 const int64_t shard_start = sharder.ShardStart(shard.Index());
341 const int64_t shard_size = sharder.ShardSize(shard.Index());
342 for (int64_t
index = shard_start;
index < shard_start + shard_size;
344 if (problem.NormWeight(
index) <= 0.0)
return false;
348 CHECK(norm_weights_are_positive);
350 if (target_radius == 0.0) {
351 return {.solution_step_size = 0.0, .objective_value = 0.0};
354 const bool objective_is_all_zeros =
355 sharder.ParallelTrueForAllShards([&](
const Sharder::Shard& shard) {
356 const int64_t shard_start = sharder.ShardStart(shard.Index());
357 const int64_t shard_size = sharder.ShardSize(shard.Index());
358 for (int64_t
index = shard_start;
index < shard_start + shard_size;
360 if (problem.Objective(
index) != 0.0)
return false;
364 if (objective_is_all_zeros) {
365 return {.solution_step_size = 0.0, .objective_value = 0.0};
368 InitialState initial_state = ComputeInitialState(problem, sharder);
372 double fixed_radius_squared = 0.0;
378 double variable_radius_coefficient =
379 initial_state.radius_coefficient_of_decided_components;
385 std::move(initial_state.undecided_components_by_shard));
392 int64_t actual_elements_seen = sharder.NumElements();
393 int64_t worst_case_elements_seen = sharder.NumElements();
396 worst_case_elements_seen +=
399 actual_elements_seen +=
402 const double step_size_threshold =
404 const double radius_squared_of_undecided_components =
406 problem, step_size_threshold, sharder,
409 const double radius_squared_at_threshold =
410 radius_squared_of_undecided_components + fixed_radius_squared +
421 VLOG(1) <<
"Total passes through variables: "
422 << actual_elements_seen /
static_cast<double>(sharder.NumElements());
423 VLOG(1) <<
"Theoretical slowdown because of shard imbalance: "
424 <<
static_cast<double>(worst_case_elements_seen) /
425 actual_elements_seen -
432 double step_size = 0.0;
433 if (variable_radius_coefficient > 0.0) {
436 variable_radius_coefficient);
443 step_size = std::numeric_limits<double>::max();
447 .solution_step_size = step_size,
454 const VectorXd& variable_lower_bounds,
455 const VectorXd& variable_upper_bounds,
456 const VectorXd& center_point,
457 const VectorXd& norm_weights,
458 const double target_radius,
460 VectorTrustRegionProblem problem(&objective_vector, &variable_lower_bounds,
461 &variable_upper_bounds, ¢er_point,
463 TrustRegionResultStepSize
solution =
464 SolveTrustRegionStepSize(problem, target_radius, sharder);
466 .solution_step_size =
solution.solution_step_size,
467 .objective_value =
solution.objective_value,
469 ComputeSolution(problem,
solution.solution_step_size, sharder),
485 const VectorXd* objective_matrix_diagonal,
488 const VectorXd* center_point,
489 const VectorXd* norm_weights)
490 : objective_vector_(*objective_vector),
491 objective_matrix_diagonal_(*objective_matrix_diagonal),
494 center_point_(*center_point),
495 norm_weight_(*norm_weights) {}
502 return variable_lower_bounds_[
index];
506 return variable_upper_bounds_[
index];
512 return objective_matrix_diagonal_[
index];
516 const VectorXd& objective_vector_;
517 const VectorXd& objective_matrix_diagonal_;
518 const VectorXd& variable_lower_bounds_;
519 const VectorXd& variable_upper_bounds_;
520 const VectorXd& center_point_;
521 const VectorXd& norm_weight_;
538class DiagonalTrustRegionProblemFromQp {
544 const VectorXd* dual_solution,
545 const VectorXd* primal_gradient,
546 const VectorXd* dual_gradient,
547 const double primal_weight)
550 dual_solution_(*dual_solution),
551 primal_gradient_(*primal_gradient),
552 dual_gradient_(*dual_gradient),
554 primal_weight_(primal_weight) {}
556 double CenterPoint(int64_t
index)
const {
557 return (
index < primal_size_) ? primal_solution_[
index]
558 : dual_solution_[
index - primal_size_];
561 double NormWeight(int64_t
index)
const {
562 return (
index < primal_size_) ? 0.5 * primal_weight_ : 0.5 / primal_weight_;
565 double LowerBound(int64_t
index)
const {
566 if (
index < primal_size_) {
567 return qp_.variable_lower_bounds[
index];
569 return std::isfinite(qp_.constraint_upper_bounds[
index - primal_size_])
570 ? -std::numeric_limits<double>::infinity()
575 double UpperBound(int64_t
index)
const {
576 if (
index < primal_size_) {
577 return qp_.variable_upper_bounds[
index];
579 return std::isfinite(qp_.constraint_lower_bounds[
index - primal_size_])
580 ? std::numeric_limits<double>::infinity()
585 double Objective(int64_t
index)
const {
586 return (
index < primal_size_) ? primal_gradient_[
index]
587 : -dual_gradient_[
index - primal_size_];
590 double ObjectiveMatrixDiagonalAt(int64_t
index)
const {
591 if (qp_.objective_matrix.has_value()) {
592 return (
index < primal_size_) ? qp_.objective_matrix->diagonal()[
index]
601 const VectorXd& primal_solution_;
602 const VectorXd& dual_solution_;
603 const VectorXd& primal_gradient_;
604 const VectorXd& dual_gradient_;
605 const int64_t primal_size_;
606 const double primal_weight_;
615template <
typename DiagonalTrustRegionProblem>
618 const double scaling_factor) {
635template <
typename DiagonalTrustRegionProblem>
638 const double scaling_factor) {
639 const double squared_norm =
642 const int64_t shard_end =
645 for (int64_t i = shard_start; i < shard_end; ++i) {
646 const double projected_coordinate =
652 return std::sqrt(squared_norm);
662template <
typename DiagonalTrustRegionProblem>
664 const Sharder& sharder,
const double target_radius,
665 const double solve_tol) {
668 double scaling_factor_lower_bound = 0.0;
669 double scaling_factor_upper_bound = 1.0;
672 scaling_factor_lower_bound = scaling_factor_upper_bound;
673 scaling_factor_upper_bound *= 2;
676 while ((scaling_factor_upper_bound - scaling_factor_lower_bound) >=
677 solve_tol * std::max(1.0, scaling_factor_lower_bound)) {
678 const double middle =
679 (scaling_factor_lower_bound + scaling_factor_upper_bound) / 2.0;
682 scaling_factor_upper_bound = middle;
684 scaling_factor_lower_bound = middle;
687 return (scaling_factor_upper_bound + scaling_factor_lower_bound) / 2.0;
693template <
typename DiagonalTrustRegionProblem>
696 const double target_radius,
const double solve_tol) {
697 CHECK_GE(target_radius, 0.0);
698 const bool norm_weights_are_positive =
701 for (int64_t i = shard_start;
709 CHECK(norm_weights_are_positive);
710 if (target_radius == 0.0) {
715 for (int64_t i = shard_start; i < shard_start + shard_size; ++i) {
719 return {.solution_step_size = 0.0,
720 .objective_value = 0.0,
723 const double optimal_scaling =
729 for (int64_t i = shard_start; i < shard_start + shard_size; ++i) {
731 const double projected_value =
737 const double final_objective_value =
739 double local_sum = 0.0;
741 for (int64_t i = shard_start;
750 return {.solution_step_size = optimal_scaling,
751 .objective_value = final_objective_value,
756 const VectorXd& objective_vector,
const VectorXd& objective_matrix_diagonal,
757 const VectorXd& variable_lower_bounds,
758 const VectorXd& variable_upper_bounds,
const VectorXd& center_point,
759 const VectorXd& norm_weights,
const double target_radius,
760 const Sharder& sharder,
const double solve_tolerance) {
762 &objective_vector, &objective_matrix_diagonal, &variable_lower_bounds,
763 &variable_upper_bounds, ¢er_point, &norm_weights);
770 const VectorXd& dual_solution,
const VectorXd& primal_gradient,
771 const VectorXd& dual_gradient,
const double primal_weight,
772 double target_radius,
const double solve_tolerance) {
775 &dual_solution, &primal_gradient,
776 &dual_gradient, primal_weight);
779 const bool norm_weights_are_positive =
782 for (int64_t i = shard_start;
784 if (problem.NormWeight(i) <= 0) {
790 CHECK(norm_weights_are_positive);
797struct MaxNormBoundResult {
812MaxNormBoundResult ComputeMaxNormPrimalTrustRegionBound(
813 const ShardedQuadraticProgram& sharded_qp,
const VectorXd&
primal_solution,
814 const double primal_radius,
const VectorXd& dual_product) {
815 LagrangianPart primal_part =
817 internal::PrimalTrustRegionProblem primal_problem(
819 TrustRegionResultStepSize trust_region_result = SolveTrustRegionStepSize(
820 primal_problem, primal_radius, sharded_qp.PrimalSharder());
821 return {.part_of_lagrangian_value = primal_part.value,
822 .trust_region_objective_delta = trust_region_result.objective_value};
825MaxNormBoundResult ComputeMaxNormDualTrustRegionBound(
826 const ShardedQuadraticProgram& sharded_qp,
const VectorXd& dual_solution,
827 const double dual_radius,
const VectorXd& primal_product) {
828 LagrangianPart dual_part =
830 internal::DualTrustRegionProblem dual_problem(
831 &sharded_qp.Qp(), &dual_solution, &dual_part.gradient);
832 TrustRegionResultStepSize trust_region_result = SolveTrustRegionStepSize(
833 dual_problem, dual_radius, sharded_qp.DualSharder());
834 return {.part_of_lagrangian_value = dual_part.value,
835 .trust_region_objective_delta = -trust_region_result.objective_value};
841double MaximumPrimalDistanceGivenWeightedDistance(
842 const double weighted_distance,
const double primal_weight) {
843 return std::sqrt(2) * weighted_distance / std::sqrt(primal_weight);
850double MaximumDualDistanceGivenWeightedDistance(
const double weighted_distance,
851 const double primal_weight) {
852 return std::sqrt(2) * weighted_distance * std::sqrt(primal_weight);
855LocalizedLagrangianBounds ComputeMaxNormLocalizedLagrangianBounds(
856 const ShardedQuadraticProgram& sharded_qp,
const VectorXd&
primal_solution,
857 const VectorXd& dual_solution,
const double primal_weight,
858 const double radius,
const Eigen::VectorXd& primal_product,
859 const Eigen::VectorXd& dual_product) {
860 const double primal_radius =
861 MaximumPrimalDistanceGivenWeightedDistance(radius, primal_weight);
862 const double dual_radius =
863 MaximumDualDistanceGivenWeightedDistance(radius, primal_weight);
868 MaxNormBoundResult primal_result = ComputeMaxNormPrimalTrustRegionBound(
871 MaxNormBoundResult dual_result = ComputeMaxNormDualTrustRegionBound(
872 sharded_qp, dual_solution, dual_radius, primal_product);
874 const double lagrangian_value = primal_result.part_of_lagrangian_value +
875 dual_result.part_of_lagrangian_value;
877 return LocalizedLagrangianBounds{
878 .lagrangian_value = lagrangian_value,
880 lagrangian_value + primal_result.trust_region_objective_delta,
882 lagrangian_value + dual_result.trust_region_objective_delta,
886LocalizedLagrangianBounds ComputeEuclideanNormLocalizedLagrangianBounds(
887 const ShardedQuadraticProgram& sharded_qp,
const VectorXd&
primal_solution,
888 const VectorXd& dual_solution,
const double primal_weight,
889 const double radius,
const Eigen::VectorXd& primal_product,
890 const Eigen::VectorXd& dual_product,
891 const bool use_diagonal_qp_trust_region_solver,
892 const double diagonal_qp_trust_region_solver_tolerance) {
893 const QuadraticProgram& qp = sharded_qp.Qp();
894 const LagrangianPart primal_part =
896 const LagrangianPart dual_part =
899 VectorXd trust_region_solution;
900 const double lagrangian_value = primal_part.value + dual_part.value;
902 Sharder joint_sharder(
903 sharded_qp.PrimalSharder(),
904 sharded_qp.PrimalSize() + sharded_qp.DualSize());
906 if (use_diagonal_qp_trust_region_solver) {
907 DiagonalTrustRegionProblemFromQp problem(
909 &dual_part.gradient, primal_weight);
912 problem, joint_sharder, radius,
913 diagonal_qp_trust_region_solver_tolerance)
916 JointTrustRegionProblem joint_problem(&qp, &
primal_solution, &dual_solution,
917 &primal_part.gradient,
918 &dual_part.gradient, primal_weight);
920 TrustRegionResultStepSize trust_region_result =
921 SolveTrustRegionStepSize(joint_problem, radius, joint_sharder);
923 trust_region_solution = ComputeSolution(
924 joint_problem, trust_region_result.solution_step_size, joint_sharder);
927 auto primal_trust_region_solution =
928 trust_region_solution.segment(0, sharded_qp.PrimalSize());
929 auto dual_trust_region_solution = trust_region_solution.segment(
930 sharded_qp.PrimalSize(), sharded_qp.DualSize());
933 double primal_objective_delta =
934 sharded_qp.PrimalSharder().ParallelSumOverShards(
935 [&](
const Sharder::Shard& shard) {
936 return shard(primal_part.gradient)
937 .dot(shard(primal_trust_region_solution) -
943 if (use_diagonal_qp_trust_region_solver &&
944 sharded_qp.Qp().objective_matrix.has_value()) {
945 primal_objective_delta += sharded_qp.PrimalSharder().ParallelSumOverShards(
946 [&](
const Sharder::Shard& shard) {
947 const int shard_start =
948 sharded_qp.PrimalSharder().ShardStart(shard.Index());
949 const int shard_size =
950 sharded_qp.PrimalSharder().ShardSize(shard.Index());
952 for (
int i = shard_start;
i < shard_start + shard_size; ++
i) {
953 sum += 0.5 * sharded_qp.Qp().objective_matrix->diagonal()[
i] *
962 const double dual_objective_delta =
963 sharded_qp.DualSharder().ParallelSumOverShards(
964 [&](
const Sharder::Shard& shard) {
965 return shard(dual_part.gradient)
966 .dot(shard(dual_trust_region_solution) - shard(dual_solution));
969 return LocalizedLagrangianBounds{
970 .lagrangian_value = lagrangian_value,
971 .lower_bound = lagrangian_value + primal_objective_delta,
972 .upper_bound = lagrangian_value + dual_objective_delta,
980 const VectorXd& dual_solution,
const PrimalDualNorm primal_dual_norm,
981 const double primal_weight,
const double radius,
982 const VectorXd* primal_product,
const VectorXd* dual_product,
983 const bool use_diagonal_qp_trust_region_solver,
984 const double diagonal_qp_trust_region_solver_tolerance) {
986 VectorXd primal_product_storage;
987 VectorXd dual_product_storage;
989 if (primal_product ==
nullptr) {
993 primal_product = &primal_product_storage;
995 if (dual_product ==
nullptr) {
996 dual_product_storage =
999 dual_product = &dual_product_storage;
1002 switch (primal_dual_norm) {
1004 return ComputeMaxNormLocalizedLagrangianBounds(
1006 *primal_product, *dual_product);
1008 return ComputeEuclideanNormLocalizedLagrangianBounds(
1010 *primal_product, *dual_product, use_diagonal_qp_trust_region_solver,
1011 diagonal_qp_trust_region_solver_tolerance);
1013 LOG(FATAL) <<
"Unrecognized primal dual norm";
static T Square(const T x)
Returns the square of x.
double UpperBound(int64_t index) const
double CenterPoint(int64_t index) const
double NormWeight(int64_t index) const
double ObjectiveMatrixDiagonalAt(int64_t index) const
double LowerBound(int64_t index) const
double Objective(int64_t index) const
DiagonalTrustRegionProblem(const VectorXd *objective_vector, const VectorXd *objective_matrix_diagonal, const VectorXd *lower_bounds, const VectorXd *upper_bounds, const VectorXd *center_point, const VectorXd *norm_weights)
int64_t PrimalSize() const
const Sharder & ConstraintMatrixSharder() const
Returns a Sharder intended for the columns of the QP's constraint matrix.
const QuadraticProgram & Qp() const
const Sharder & TransposedConstraintMatrixSharder() const
Returns a Sharder intended for the rows of the QP's constraint matrix.
const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > & TransposedConstraintMatrix() const
Returns a reference to the transpose of the QP's constraint matrix.
const Sharder & PrimalSharder() const
Returns a Sharder intended for primal vectors.
bool ParallelTrueForAllShards(const std::function< bool(const Shard &)> &func) const
double ParallelSumOverShards(const std::function< double(const Shard &)> &func) const
Runs func on each of the shards and sums the results.
int64_t ShardSize(int shard) const
void ParallelForEachShard(const std::function< void(const Shard &)> &func) const
Runs func on each of the shards.
int64_t NumElements() const
The number of elements that are split into shards.
int64_t ShardStart(int shard) const
double LowerBound(const LinearExpression &linear_expression)
double UpperBound(const LinearExpression &linear_expression)
double ComputeInitialUndecidedComponents(const TrustRegionProblem &problem, int64_t start_index, int64_t end_index, std::vector< int64_t > &undecided_components)
double ProjectedValue(const TrustRegionProblem &problem, const int64_t index, const double step_size)
double RemoveCriticalStepsAboveThreshold(const TrustRegionProblem &problem, const double step_size_threshold, std::vector< int64_t > &undecided_components)
double CriticalStepSize(const TrustRegionProblem &problem, const int64_t index)
double RemoveCriticalStepsBelowThreshold(const TrustRegionProblem &problem, const double step_size_threshold, std::vector< int64_t > &undecided_components)
double EasyMedian(ArrayType array, ValueFunction value_function)
double RadiusSquaredOfUndecidedComponents(const TrustRegionProblem &problem, const double step_size, const std::vector< int64_t > &undecided_components)
Validation utilities for solvers.proto.
LagrangianPart ComputeDualGradient(const ShardedQuadraticProgram &sharded_qp, const VectorXd &dual_solution, const VectorXd &primal_product)
LagrangianPart ComputePrimalGradient(const ShardedQuadraticProgram &sharded_qp, const VectorXd &primal_solution, const VectorXd &dual_product)
TrustRegionResult SolveDiagonalTrustRegionProblem(const DiagonalTrustRegionProblem &problem, const Sharder &sharder, const double target_radius, const double solve_tol)
double NormOfDeltaProjection(const DiagonalTrustRegionProblem &problem, const Sharder &sharder, const double scaling_factor)
VectorXd TransposedMatrixVectorProduct(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &vector, const Sharder &sharder)
TrustRegionResult SolveDiagonalTrustRegion(const VectorXd &objective_vector, const VectorXd &objective_matrix_diagonal, const VectorXd &variable_lower_bounds, const VectorXd &variable_upper_bounds, const VectorXd ¢er_point, const VectorXd &norm_weights, const double target_radius, const Sharder &sharder, const double solve_tolerance)
TrustRegionResult SolveDiagonalQpTrustRegion(const ShardedQuadraticProgram &sharded_qp, const VectorXd &primal_solution, const VectorXd &dual_solution, const VectorXd &primal_gradient, const VectorXd &dual_gradient, const double primal_weight, double target_radius, const double solve_tolerance)
double FindScalingFactor(const DiagonalTrustRegionProblem &problem, const Sharder &sharder, const double target_radius, const double solve_tol)
@ kMaxNorm
The joint norm ||(x,y)||_PD = max{|x|_P, |y|_D}.
@ kEuclideanNorm
The joint norm (||(x,y)||_PD)^2 = (|x|_P)^2 + (|y|_D)^2.
TrustRegionResult SolveTrustRegion(const VectorXd &objective_vector, const VectorXd &variable_lower_bounds, const VectorXd &variable_upper_bounds, const VectorXd ¢er_point, const VectorXd &norm_weights, const double target_radius, const Sharder &sharder)
double ProjectedValueOfScaledDifference(const DiagonalTrustRegionProblem &problem, const int64_t index, const double scaling_factor)
LocalizedLagrangianBounds ComputeLocalizedLagrangianBounds(const ShardedQuadraticProgram &sharded_qp, const VectorXd &primal_solution, const VectorXd &dual_solution, const PrimalDualNorm primal_dual_norm, const double primal_weight, const double radius, const VectorXd *primal_product, const VectorXd *dual_product, const bool use_diagonal_qp_trust_region_solver, const double diagonal_qp_trust_region_solver_tolerance)
Coefficient ComputeObjectiveValue(const LinearBooleanProblem &problem, const std::vector< bool > &assignment)
Returns the objective value under the current assignment.
std::vector< double > lower_bounds
std::vector< double > upper_bounds
Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > constraint_matrix
Eigen::VectorXd solution
The solution.
double trust_region_objective_delta
double part_of_lagrangian_value
double objective_value
The value objective_vector^T * (solution - center_point).
std::vector< std::vector< int64_t > > undecided_components_by_shard
double solution_step_size
The step_size of the solution.
double radius_coefficient_of_decided_components