Google OR-Tools v9.12
a fast and portable software suite for combinatorial optimization
Loading...
Searching...
No Matches
sharded_quadratic_program.h
Go to the documentation of this file.
1// Copyright 2010-2025 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
14#ifndef PDLP_SHARDED_QUADRATIC_PROGRAM_H_
15#define PDLP_SHARDED_QUADRATIC_PROGRAM_H_
16
17#include <cstdint>
18#include <memory>
19#include <optional>
20
21#include "Eigen/Core"
22#include "Eigen/SparseCore"
26#include "ortools/pdlp/solvers.pb.h"
28
30
31// This class stores:
32// - A `QuadraticProgram` (QP)
33// - A transposed version of the QP's constraint matrix
34// - A thread scheduler
35// - Various `Sharder` objects for doing sharded matrix and vector
36// computations.
38 public:
39 // Requires `num_shards` >= `num_threads` >= 1.
40 // Note that the `qp` is intentionally passed by value.
41 // If `logger` is not nullptr, warns about unbalanced matrices using it;
42 // otherwise warns via Google standard logging.
44 QuadraticProgram qp, int num_threads, int num_shards,
45 SchedulerType scheduler_type = SCHEDULER_TYPE_GOOGLE_THREADPOOL,
46 operations_research::SolverLogger* logger = nullptr);
47
48 // Movable but not copyable.
53
54 const QuadraticProgram& Qp() const { return qp_; }
55
56 // Returns a reference to the transpose of the QP's constraint matrix.
57 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>&
59 return transposed_constraint_matrix_;
60 }
61
62 // Returns a `Sharder` intended for the columns of the QP's constraint matrix.
64 return constraint_matrix_sharder_;
65 }
66 // Returns a `Sharder` intended for the rows of the QP's constraint matrix.
68 return transposed_constraint_matrix_sharder_;
69 }
70 // Returns a `Sharder` intended for primal vectors.
71 const Sharder& PrimalSharder() const { return primal_sharder_; }
72 // Returns a `Sharder` intended for dual vectors.
73 const Sharder& DualSharder() const { return dual_sharder_; }
74
75 int64_t PrimalSize() const { return qp_.variable_lower_bounds.size(); }
76 int64_t DualSize() const { return qp_.constraint_lower_bounds.size(); }
77
78 // Rescales the QP (including objective, variable bounds, constraint bounds,
79 // constraint matrix, and transposed constraint matrix) based on
80 // `col_scaling_vec` and `row_scaling_vec`. That is, rescale the problem so
81 // that each variable is rescaled as variable[i] <- variable[i] /
82 // `col_scaling_vec[i]`, and the j-th constraint is multiplied by
83 // `row_scaling_vec[j]`. `col_scaling_vec` and `row_scaling_vec` must be
84 // positive.
85 void RescaleQuadraticProgram(const Eigen::VectorXd& col_scaling_vec,
86 const Eigen::VectorXd& row_scaling_vec);
87
88 void SwapVariableBounds(Eigen::VectorXd& variable_lower_bounds,
89 Eigen::VectorXd& variable_upper_bounds) {
90 qp_.variable_lower_bounds.swap(variable_lower_bounds);
91 qp_.variable_upper_bounds.swap(variable_upper_bounds);
92 }
93
94 void SwapConstraintBounds(Eigen::VectorXd& constraint_lower_bounds,
95 Eigen::VectorXd& constraint_upper_bounds) {
96 qp_.constraint_lower_bounds.swap(constraint_lower_bounds);
97 qp_.constraint_upper_bounds.swap(constraint_upper_bounds);
98 }
99
100 void SetConstraintBounds(int64_t constraint_index,
101 std::optional<double> lower_bound,
102 std::optional<double> upper_bound);
103
104 // Swaps `objective` with the `objective_vector` in the quadratic program.
105 // Swapping `objective_matrix` is not yet supported because it hasn't been
106 // needed.
107 void SwapObjectiveVector(Eigen::VectorXd& objective) {
108 qp_.objective_vector.swap(objective);
109 }
110
111 // Replaces constraint bounds whose absolute value is >= `threshold` with
112 // the corresponding infinity.
113 void ReplaceLargeConstraintBoundsWithInfinity(double threshold);
114
115 private:
117 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>
118 transposed_constraint_matrix_;
119 std::unique_ptr<Scheduler> scheduler_;
120 Sharder constraint_matrix_sharder_;
121 Sharder transposed_constraint_matrix_sharder_;
122 Sharder primal_sharder_;
123 Sharder dual_sharder_;
124};
125
126} // namespace operations_research::pdlp
127
128#endif // PDLP_SHARDED_QUADRATIC_PROGRAM_H_
ShardedQuadraticProgram(ShardedQuadraticProgram &&)=default
const Sharder & ConstraintMatrixSharder() const
Returns a Sharder intended for the columns of the QP's constraint matrix.
void SwapConstraintBounds(Eigen::VectorXd &constraint_lower_bounds, Eigen::VectorXd &constraint_upper_bounds)
const Sharder & TransposedConstraintMatrixSharder() const
Returns a Sharder intended for the rows of the QP's constraint matrix.
void SwapVariableBounds(Eigen::VectorXd &variable_lower_bounds, Eigen::VectorXd &variable_upper_bounds)
const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > & TransposedConstraintMatrix() const
Returns a reference to the transpose of the QP's constraint matrix.
ShardedQuadraticProgram(const ShardedQuadraticProgram &)=delete
Movable but not copyable.
void SetConstraintBounds(int64_t constraint_index, std::optional< double > lower_bound, std::optional< double > upper_bound)
ShardedQuadraticProgram & operator=(ShardedQuadraticProgram &&)=default
const Sharder & DualSharder() const
Returns a Sharder intended for dual vectors.
ShardedQuadraticProgram & operator=(const ShardedQuadraticProgram &)=delete
ShardedQuadraticProgram(QuadraticProgram qp, int num_threads, int num_shards, SchedulerType scheduler_type=SCHEDULER_TYPE_GOOGLE_THREADPOOL, operations_research::SolverLogger *logger=nullptr)
const Sharder & PrimalSharder() const
Returns a Sharder intended for primal vectors.
void RescaleQuadraticProgram(const Eigen::VectorXd &col_scaling_vec, const Eigen::VectorXd &row_scaling_vec)
Validation utilities for solvers.proto.