24#include "absl/types/optional.h"
34using ::Eigen::VectorXd;
74class VectorTrustRegionProblem {
76 VectorTrustRegionProblem(
const VectorXd* objective,
80 const VectorXd* norm_weight)
85 norm_weight_(*norm_weight) {}
89 double CenterPoint(int64_t
index)
const {
return center_point_(
index); }
90 double NormWeight(int64_t
index)
const {
return norm_weight_(
index); }
94 const VectorXd& lower_bound_;
95 const VectorXd& upper_bound_;
96 const VectorXd& center_point_;
97 const VectorXd& norm_weight_;
114class JointTrustRegionProblem {
116 JointTrustRegionProblem(
const QuadraticProgram* qp,
117 const VectorXd* primal_solution,
118 const VectorXd* dual_solution,
119 const VectorXd* primal_gradient,
120 const VectorXd* dual_gradient,
121 const double primal_weight)
124 primal_solution_(*primal_solution),
125 dual_solution_(*dual_solution),
126 primal_gradient_(*primal_gradient),
127 dual_gradient_(*dual_gradient),
128 primal_weight_(primal_weight) {}
129 double Objective(int64_t
index)
const {
130 return index < primal_size_ ? primal_gradient_[
index]
131 : -dual_gradient_[
index - primal_size_];
134 return index < primal_size_ ? qp_.variable_lower_bounds[
index]
135 : std::isfinite(qp_.constraint_upper_bounds[
index - primal_size_])
136 ? -std::numeric_limits<double>::infinity()
140 return index < primal_size_ ? qp_.variable_upper_bounds[
index]
141 : std::isfinite(qp_.constraint_lower_bounds[
index - primal_size_])
142 ? std::numeric_limits<double>::infinity()
145 double CenterPoint(int64_t
index)
const {
146 return index < primal_size_ ? primal_solution_[
index]
147 : dual_solution_[
index - primal_size_];
149 double NormWeight(int64_t
index)
const {
150 return index < primal_size_ ? 0.5 * primal_weight_ : 0.5 / primal_weight_;
154 const QuadraticProgram& qp_;
155 const int64_t primal_size_;
156 const VectorXd& primal_solution_;
157 const VectorXd& dual_solution_;
158 const VectorXd& primal_gradient_;
159 const VectorXd& dual_gradient_;
160 const double primal_weight_;
163struct TrustRegionResultStepSize {
176template <
typename TrustRegionProblem>
177double MedianOfShardMedians(
178 const TrustRegionProblem& problem,
179 const std::vector<std::vector<int64_t>>& indexed_components_by_shard,
180 const Sharder& sharder) {
181 std::vector<absl::optional<double>> shard_medians(sharder.NumShards(),
183 sharder.ParallelForEachShard([&](
const Sharder::Shard& shard) {
184 const auto& indexed_shard_components =
185 indexed_components_by_shard[shard.Index()];
186 if (!indexed_shard_components.empty()) {
187 shard_medians[shard.Index()] = internal::EasyMedian(
188 indexed_shard_components, [&](const int64_t index) {
189 return internal::CriticalStepSize(problem, index);
193 std::vector<double> non_empty_medians;
194 for (
const auto& median : shard_medians) {
195 if (median.has_value()) {
196 non_empty_medians.push_back(*median);
199 CHECK(!non_empty_medians.empty());
201 [](
const double x) {
return x; });
209template <
typename TrustRegionProblem>
210InitialState ComputeInitialState(
const TrustRegionProblem& problem,
211 const Sharder& sharder) {
213 result.undecided_components_by_shard.resize(sharder.NumShards());
214 result.radius_coefficient_of_decided_components =
215 sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
216 const int64_t shard_start = sharder.ShardStart(shard.Index());
217 const int64_t shard_size = sharder.ShardSize(shard.Index());
219 problem, shard_start, shard_start + shard_size,
220 result.undecided_components_by_shard[shard.Index()]);
225template <
typename TrustRegionProblem>
227 const TrustRegionProblem& problem,
const double step_size,
228 const Sharder& sharder,
230 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
236template <
typename TrustRegionProblem>
238 const TrustRegionProblem& problem,
const double step_size_threshold,
239 const Sharder& sharder,
241 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
243 problem, step_size_threshold,
248template <
typename TrustRegionProblem>
250 const TrustRegionProblem& problem,
const double step_size_threshold,
251 const Sharder& sharder,
253 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
255 problem, step_size_threshold,
260int64_t NumUndecidedComponents(
262 int64_t num_undecided_components = 0;
264 num_undecided_components += undecided_components.size();
266 return num_undecided_components;
269int64_t MaxUndecidedComponentsInAnyShard(
273 max = std::max<int64_t>(
max, undecided_components.size());
278template <
typename TrustRegionProblem>
279VectorXd ComputeSolution(
const TrustRegionProblem& problem,
280 const double step_size,
const Sharder& sharder) {
281 VectorXd solution(sharder.NumElements());
282 sharder.ParallelForEachShard([&](
const Sharder::Shard& shard) {
283 const int64_t shard_start = sharder.ShardStart(shard.Index());
284 const int64_t shard_size = sharder.ShardSize(shard.Index());
285 for (int64_t
index = shard_start;
index < shard_start + shard_size;
293template <
typename TrustRegionProblem>
295 const double step_size,
const Sharder& sharder) {
296 return sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
297 const int64_t shard_start = sharder.ShardStart(shard.Index());
298 const int64_t shard_size = sharder.ShardSize(shard.Index());
299 double shard_value = 0.0;
300 for (int64_t
index = shard_start;
index < shard_start + shard_size;
302 shard_value += problem.Objective(
index) *
304 problem.CenterPoint(
index));
331template <
typename TrustRegionProblem>
332TrustRegionResultStepSize SolveTrustRegionStepSize(
333 const TrustRegionProblem& problem,
const double target_radius,
334 const Sharder& sharder) {
337 const bool norm_weights_are_positive =
338 sharder.ParallelTrueForAllShards([&](
const Sharder::Shard& shard) {
339 const int64_t shard_start = sharder.ShardStart(shard.Index());
340 const int64_t shard_size = sharder.ShardSize(shard.Index());
341 for (int64_t
index = shard_start;
index < shard_start + shard_size;
343 if (problem.NormWeight(
index) <= 0.0)
return false;
347 CHECK(norm_weights_are_positive);
349 if (target_radius == 0.0) {
350 return {.solution_step_size = 0.0, .objective_value = 0.0};
353 const bool objective_is_all_zeros =
354 sharder.ParallelTrueForAllShards([&](
const Sharder::Shard& shard) {
355 const int64_t shard_start = sharder.ShardStart(shard.Index());
356 const int64_t shard_size = sharder.ShardSize(shard.Index());
357 for (int64_t
index = shard_start;
index < shard_start + shard_size;
359 if (problem.Objective(
index) != 0.0)
return false;
363 if (objective_is_all_zeros) {
364 return {.solution_step_size = 0.0, .objective_value = 0.0};
367 InitialState initial_state = ComputeInitialState(problem, sharder);
371 double fixed_radius_squared = 0.0;
377 double variable_radius_coefficient =
378 initial_state.radius_coefficient_of_decided_components;
384 std::move(initial_state.undecided_components_by_shard));
391 int64_t actual_elements_seen = sharder.NumElements();
392 int64_t worst_case_elements_seen = sharder.NumElements();
395 worst_case_elements_seen +=
398 actual_elements_seen +=
401 const double step_size_threshold =
403 const double radius_squared_of_undecided_components =
405 problem, step_size_threshold, sharder,
408 const double radius_squared_at_threshold =
409 radius_squared_of_undecided_components + fixed_radius_squared +
420 VLOG(1) <<
"Total passes through variables: "
421 << actual_elements_seen /
static_cast<double>(sharder.NumElements());
422 VLOG(1) <<
"Theoretical slowdown because of shard imbalance: "
423 <<
static_cast<double>(worst_case_elements_seen) /
424 actual_elements_seen -
431 double step_size = 0.0;
432 if (variable_radius_coefficient > 0.0) {
435 variable_radius_coefficient);
446 .solution_step_size = step_size,
457 const double target_radius,
462 TrustRegionResultStepSize solution =
463 SolveTrustRegionStepSize(problem, target_radius, sharder);
466 .objective_value = solution.objective_value,
468 ComputeSolution(problem, solution.solution_step_size, sharder),
501 return variable_lower_bounds_[
index];
505 return variable_upper_bounds_[
index];
511 return objective_matrix_diagonal_[
index];
515 const VectorXd& objective_vector_;
516 const VectorXd& objective_matrix_diagonal_;
517 const VectorXd& variable_lower_bounds_;
518 const VectorXd& variable_upper_bounds_;
519 const VectorXd& center_point_;
520 const VectorXd& norm_weight_;
542 const VectorXd* primal_solution,
543 const VectorXd* dual_solution,
544 const VectorXd* primal_gradient,
545 const VectorXd* dual_gradient,
546 const double primal_weight)
548 primal_solution_(*primal_solution),
549 dual_solution_(*dual_solution),
550 primal_gradient_(*primal_gradient),
551 dual_gradient_(*dual_gradient),
552 primal_size_(primal_solution->size()),
553 primal_weight_(primal_weight) {}
556 return (
index < primal_size_) ? primal_solution_[
index]
557 : dual_solution_[
index - primal_size_];
561 return (
index < primal_size_) ? 0.5 * primal_weight_ : 0.5 / primal_weight_;
565 if (
index < primal_size_) {
566 return qp_.variable_lower_bounds[
index];
568 return std::isfinite(qp_.constraint_upper_bounds[
index - primal_size_])
569 ? -std::numeric_limits<double>::infinity()
575 if (
index < primal_size_) {
576 return qp_.variable_upper_bounds[
index];
578 return std::isfinite(qp_.constraint_lower_bounds[
index - primal_size_])
579 ? std::numeric_limits<double>::infinity()
585 return (
index < primal_size_) ? primal_gradient_[
index]
586 : -dual_gradient_[
index - primal_size_];
590 if (qp_.objective_matrix.has_value()) {
591 return (
index < primal_size_) ? qp_.objective_matrix->diagonal()[
index]
600 const VectorXd& primal_solution_;
601 const VectorXd& dual_solution_;
602 const VectorXd& primal_gradient_;
603 const VectorXd& dual_gradient_;
604 const int64_t primal_size_;
605 const double primal_weight_;
614template <
typename DiagonalTrustRegionProblem>
617 const double scaling_factor) {
634template <
typename DiagonalTrustRegionProblem>
637 const double scaling_factor) {
638 const double squared_norm =
641 const int64_t shard_end =
644 for (int64_t i = shard_start; i < shard_end; ++i) {
645 const double projected_coordinate =
651 return std::sqrt(squared_norm);
661template <
typename DiagonalTrustRegionProblem>
663 const Sharder& sharder,
const double target_radius,
664 const double solve_tol) {
667 double scaling_factor_lower_bound = 0.0;
668 double scaling_factor_upper_bound = 1.0;
671 scaling_factor_lower_bound = scaling_factor_upper_bound;
672 scaling_factor_upper_bound *= 2;
675 while ((scaling_factor_upper_bound - scaling_factor_lower_bound) >=
676 solve_tol *
std::max(1.0, scaling_factor_lower_bound)) {
677 const double middle =
678 (scaling_factor_lower_bound + scaling_factor_upper_bound) / 2.0;
681 scaling_factor_upper_bound = middle;
683 scaling_factor_lower_bound = middle;
686 return (scaling_factor_upper_bound + scaling_factor_lower_bound) / 2.0;
692template <
typename DiagonalTrustRegionProblem>
695 const double target_radius,
const double solve_tol) {
697 const bool norm_weights_are_positive =
700 for (int64_t i = shard_start;
702 if (problem.NormWeight(i) <= 0) {
708 CHECK(norm_weights_are_positive);
709 const double optimal_scaling =
711 VectorXd solution(sharder.NumElements());
712 sharder.ParallelForEachShard([&](
const Sharder::Shard& shard) {
713 const int64_t shard_start = sharder.ShardStart(shard.Index());
714 const int64_t shard_size = sharder.ShardSize(shard.Index());
715 for (int64_t i = shard_start; i < shard_start + shard_size; ++i) {
716 const double weight = problem.NormWeight(i);
717 const double projected_value =
720 problem.CenterPoint(i) + std::sqrt(1 /
weight) * projected_value;
723 const double final_objective_value =
724 sharder.ParallelSumOverShards([&](
const Sharder::Shard& shard) {
725 double local_sum = 0.0;
726 const int64_t shard_start = sharder.ShardStart(shard.Index());
727 for (int64_t i = shard_start;
728 i < shard_start + sharder.ShardSize(shard.Index()); ++i) {
729 const double diff = solution[i] - problem.CenterPoint(i);
731 0.5 * diff * problem.ObjectiveMatrixDiagonalAt(i) * diff +
732 diff * problem.Objective(i);
736 return {.solution_step_size = optimal_scaling,
737 .objective_value = final_objective_value,
738 .solution = solution};
745 const VectorXd&
norm_weights,
const double target_radius,
746 const Sharder& sharder,
const double solve_tolerance) {
756 const VectorXd& dual_solution,
const VectorXd& primal_gradient,
757 const VectorXd& dual_gradient,
const double primal_weight,
758 double target_radius,
const double solve_tolerance) {
761 &dual_solution, &primal_gradient,
762 &dual_gradient, primal_weight);
765 const bool norm_weights_are_positive =
768 for (int64_t i = shard_start;
770 if (problem.NormWeight(i) <= 0) {
776 CHECK(norm_weights_are_positive);
783struct MaxNormBoundResult {
798MaxNormBoundResult ComputeMaxNormPrimalTrustRegionBound(
799 const ShardedQuadraticProgram& sharded_qp,
const VectorXd& primal_solution,
800 const double primal_radius,
const VectorXd& dual_product) {
801 LagrangianPart primal_part =
803 internal::PrimalTrustRegionProblem primal_problem(
804 &sharded_qp.Qp(), &primal_solution, &primal_part.gradient);
805 TrustRegionResultStepSize trust_region_result = SolveTrustRegionStepSize(
806 primal_problem, primal_radius, sharded_qp.PrimalSharder());
807 return {.part_of_lagrangian_value = primal_part.value,
808 .trust_region_objective_delta = trust_region_result.objective_value};
811MaxNormBoundResult ComputeMaxNormDualTrustRegionBound(
812 const ShardedQuadraticProgram& sharded_qp,
const VectorXd& dual_solution,
813 const double dual_radius,
const VectorXd& primal_product) {
814 LagrangianPart dual_part =
816 internal::DualTrustRegionProblem dual_problem(
817 &sharded_qp.Qp(), &dual_solution, &dual_part.gradient);
818 TrustRegionResultStepSize trust_region_result = SolveTrustRegionStepSize(
819 dual_problem, dual_radius, sharded_qp.DualSharder());
820 return {.part_of_lagrangian_value = dual_part.value,
821 .trust_region_objective_delta = -trust_region_result.objective_value};
827double MaximumPrimalDistanceGivenWeightedDistance(
828 const double weighted_distance,
const double primal_weight) {
829 return std::sqrt(2) * weighted_distance / std::sqrt(primal_weight);
835double MaximumDualDistanceGivenWeightedDistance(
const double weighted_distance,
836 const double primal_weight) {
837 return std::sqrt(2) * weighted_distance * std::sqrt(primal_weight);
840LocalizedLagrangianBounds ComputeMaxNormLocalizedLagrangianBounds(
841 const ShardedQuadraticProgram& sharded_qp,
const VectorXd& primal_solution,
842 const VectorXd& dual_solution,
const double primal_weight,
843 const double radius,
const Eigen::VectorXd& primal_product,
844 const Eigen::VectorXd& dual_product) {
845 const double primal_radius =
846 MaximumPrimalDistanceGivenWeightedDistance(radius, primal_weight);
847 const double dual_radius =
848 MaximumDualDistanceGivenWeightedDistance(radius, primal_weight);
853 MaxNormBoundResult primal_result = ComputeMaxNormPrimalTrustRegionBound(
854 sharded_qp, primal_solution, primal_radius, dual_product);
856 MaxNormBoundResult dual_result = ComputeMaxNormDualTrustRegionBound(
857 sharded_qp, dual_solution, dual_radius, primal_product);
859 const double lagrangian_value = primal_result.part_of_lagrangian_value +
860 dual_result.part_of_lagrangian_value;
862 return LocalizedLagrangianBounds{
863 .lagrangian_value = lagrangian_value,
865 lagrangian_value + primal_result.trust_region_objective_delta,
867 lagrangian_value + dual_result.trust_region_objective_delta,
871LocalizedLagrangianBounds ComputeEuclideanNormLocalizedLagrangianBounds(
872 const ShardedQuadraticProgram& sharded_qp,
const VectorXd& primal_solution,
873 const VectorXd& dual_solution,
const double primal_weight,
874 const double radius,
const Eigen::VectorXd& primal_product,
875 const Eigen::VectorXd& dual_product,
876 const bool use_diagonal_qp_trust_region_solver,
877 const double diagonal_qp_trust_region_solver_tolerance) {
878 const QuadraticProgram& qp = sharded_qp.Qp();
879 const LagrangianPart primal_part =
881 const LagrangianPart dual_part =
884 VectorXd trust_region_solution;
885 const double lagrangian_value = primal_part.value + dual_part.value;
887 Sharder joint_sharder(
888 sharded_qp.PrimalSharder(),
889 sharded_qp.PrimalSize() + sharded_qp.DualSize());
891 if (use_diagonal_qp_trust_region_solver) {
892 DiagonalTrustRegionProblemFromQp problem(
893 &qp, &primal_solution, &dual_solution, &primal_part.gradient,
894 &dual_part.gradient, primal_weight);
897 problem, joint_sharder, radius,
898 diagonal_qp_trust_region_solver_tolerance)
901 JointTrustRegionProblem joint_problem(&qp, &primal_solution, &dual_solution,
902 &primal_part.gradient,
903 &dual_part.gradient, primal_weight);
905 TrustRegionResultStepSize trust_region_result =
906 SolveTrustRegionStepSize(joint_problem, radius, joint_sharder);
908 trust_region_solution = ComputeSolution(
909 joint_problem, trust_region_result.solution_step_size, joint_sharder);
912 auto primal_trust_region_solution =
913 trust_region_solution.segment(0, sharded_qp.PrimalSize());
914 auto dual_trust_region_solution = trust_region_solution.segment(
915 sharded_qp.PrimalSize(), sharded_qp.DualSize());
918 double primal_objective_delta =
919 sharded_qp.PrimalSharder().ParallelSumOverShards(
920 [&](
const Sharder::Shard& shard) {
921 return shard(primal_part.gradient)
922 .dot(shard(primal_trust_region_solution) -
923 shard(primal_solution));
928 if (use_diagonal_qp_trust_region_solver &&
929 sharded_qp.Qp().objective_matrix.has_value()) {
930 primal_objective_delta += sharded_qp.PrimalSharder().ParallelSumOverShards(
931 [&](
const Sharder::Shard& shard) {
932 const int shard_start =
933 sharded_qp.PrimalSharder().ShardStart(shard.Index());
934 const int shard_size =
935 sharded_qp.PrimalSharder().ShardSize(shard.Index());
937 for (
int i = shard_start; i < shard_start + shard_size; ++i) {
938 sum += 0.5 * sharded_qp.Qp().objective_matrix->diagonal()[i] *
947 const double dual_objective_delta =
948 sharded_qp.DualSharder().ParallelSumOverShards(
949 [&](
const Sharder::Shard& shard) {
950 return shard(dual_part.gradient)
951 .dot(shard(dual_trust_region_solution) - shard(dual_solution));
954 return LocalizedLagrangianBounds{
955 .lagrangian_value = lagrangian_value,
956 .lower_bound = lagrangian_value + primal_objective_delta,
957 .upper_bound = lagrangian_value + dual_objective_delta,
965 const VectorXd& dual_solution,
const PrimalDualNorm primal_dual_norm,
966 const double primal_weight,
const double radius,
967 const VectorXd* primal_product,
const VectorXd* dual_product,
968 const bool use_diagonal_qp_trust_region_solver,
969 const double diagonal_qp_trust_region_solver_tolerance) {
971 VectorXd primal_product_storage;
972 VectorXd dual_product_storage;
974 if (primal_product ==
nullptr) {
978 primal_product = &primal_product_storage;
980 if (dual_product ==
nullptr) {
981 dual_product_storage =
984 dual_product = &dual_product_storage;
987 switch (primal_dual_norm) {
988 case PrimalDualNorm::kMaxNorm:
989 return ComputeMaxNormLocalizedLagrangianBounds(
990 sharded_qp, primal_solution, dual_solution, primal_weight, radius,
991 *primal_product, *dual_product);
992 case PrimalDualNorm::kEuclideanNorm:
993 return ComputeEuclideanNormLocalizedLagrangianBounds(
994 sharded_qp, primal_solution, dual_solution, primal_weight, radius,
995 *primal_product, *dual_product, use_diagonal_qp_trust_region_solver,
996 diagonal_qp_trust_region_solver_tolerance);
998 LOG(
FATAL) <<
"Unrecognized primal dual norm";
#define CHECK_GE(val1, val2)
#define VLOG(verboselevel)
static T Square(const T x)
double Objective(int64_t index) const
double UpperBound(int64_t index) const
double ObjectiveMatrixDiagonalAt(int64_t index) const
double LowerBound(int64_t index) const
DiagonalTrustRegionProblemFromQp(const QuadraticProgram *qp, const VectorXd *primal_solution, const VectorXd *dual_solution, const VectorXd *primal_gradient, const VectorXd *dual_gradient, const double primal_weight)
double CenterPoint(int64_t index) const
double NormWeight(int64_t index) const
double Objective(int64_t index) const
double UpperBound(int64_t index) const
double ObjectiveMatrixDiagonalAt(int64_t index) const
double LowerBound(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)
double CenterPoint(int64_t index) const
double NormWeight(int64_t index) const
const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > & TransposedConstraintMatrix() const
const Sharder & TransposedConstraintMatrixSharder() const
const Sharder & PrimalSharder() const
int64_t PrimalSize() const
const QuadraticProgram & Qp() const
const Sharder & ConstraintMatrixSharder() const
double ParallelSumOverShards(const std::function< double(const Shard &)> &func) const
bool ParallelTrueForAllShards(const std::function< bool(const Shard &)> &func) const
int64_t ShardSize(int shard) const
int64_t ShardStart(int shard) const
Fractional Square(Fractional f)
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 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)
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)
LagrangianPart ComputeDualGradient(const ShardedQuadraticProgram &sharded_qp, const Eigen::VectorXd &dual_solution, const Eigen::VectorXd &primal_product)
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)
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)
std::function< int64_t(const Model &)> LowerBound(IntegerVariable v)
std::function< int64_t(const Model &)> UpperBound(IntegerVariable v)
Coefficient ComputeObjectiveValue(const LinearBooleanProblem &problem, const std::vector< bool > &assignment)
std::vector< double > lower_bounds
std::vector< double > upper_bounds
Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > constraint_matrix
double solution_step_size
double trust_region_objective_delta
double part_of_lagrangian_value
std::vector< std::vector< int64_t > > undecided_components_by_shard
double solution_step_size
double radius_coefficient_of_decided_components
VectorXd variable_lower_bounds
VectorXd variable_upper_bounds
VectorXd objective_vector
VectorXd objective_matrix_diagonal