Google OR-Tools v9.11
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-2024 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#include <utility>
21
22#include "Eigen/Core"
23#include "Eigen/SparseCore"
28
30
31// This class stores:
32// - A `QuadraticProgram` (QP)
33// - A transposed version of the QP's constraint matrix
34// - A thread pool
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.
43 ShardedQuadraticProgram(QuadraticProgram qp, int num_threads, int num_shards,
44 operations_research::SolverLogger* logger = nullptr);
45
46 // Movable but not copyable.
51
52 const QuadraticProgram& Qp() const { return qp_; }
53
54 // Returns a reference to the transpose of the QP's constraint matrix.
55 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>&
57 return transposed_constraint_matrix_;
58 }
59
60 // Returns a `Sharder` intended for the columns of the QP's constraint matrix.
62 return constraint_matrix_sharder_;
63 }
64 // Returns a `Sharder` intended for the rows of the QP's constraint matrix.
66 return transposed_constraint_matrix_sharder_;
67 }
68 // Returns a `Sharder` intended for primal vectors.
69 const Sharder& PrimalSharder() const { return primal_sharder_; }
70 // Returns a `Sharder` intended for dual vectors.
71 const Sharder& DualSharder() const { return dual_sharder_; }
72
73 int64_t PrimalSize() const { return qp_.variable_lower_bounds.size(); }
74 int64_t DualSize() const { return qp_.constraint_lower_bounds.size(); }
75
76 // Rescales the QP (including objective, variable bounds, constraint bounds,
77 // constraint matrix, and transposed constraint matrix) based on
78 // `col_scaling_vec` and `row_scaling_vec`. That is, rescale the problem so
79 // that each variable is rescaled as variable[i] <- variable[i] /
80 // `col_scaling_vec[i]`, and the j-th constraint is multiplied by
81 // `row_scaling_vec[j]`. `col_scaling_vec` and `row_scaling_vec` must be
82 // positive.
83 void RescaleQuadraticProgram(const Eigen::VectorXd& col_scaling_vec,
84 const Eigen::VectorXd& row_scaling_vec);
85
86 void SwapVariableBounds(Eigen::VectorXd& variable_lower_bounds,
87 Eigen::VectorXd& variable_upper_bounds) {
88 qp_.variable_lower_bounds.swap(variable_lower_bounds);
89 qp_.variable_upper_bounds.swap(variable_upper_bounds);
90 }
91
92 void SwapConstraintBounds(Eigen::VectorXd& constraint_lower_bounds,
93 Eigen::VectorXd& constraint_upper_bounds) {
94 qp_.constraint_lower_bounds.swap(constraint_lower_bounds);
95 qp_.constraint_upper_bounds.swap(constraint_upper_bounds);
96 }
97
99 std::optional<double> lower_bound,
100 std::optional<double> upper_bound);
101
102 // Swaps `objective` with the `objective_vector` in the quadratic program.
103 // Swapping `objective_matrix` is not yet supported because it hasn't been
104 // needed.
105 void SwapObjectiveVector(Eigen::VectorXd& objective) {
106 qp_.objective_vector.swap(objective);
107 }
108
109 // Replaces constraint bounds whose absolute value is >= `threshold` with
110 // the corresponding infinity.
111 void ReplaceLargeConstraintBoundsWithInfinity(double threshold);
112
113 private:
115 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>
116 transposed_constraint_matrix_;
117 std::unique_ptr<ThreadPool> thread_pool_;
118 Sharder constraint_matrix_sharder_;
119 Sharder transposed_constraint_matrix_sharder_;
120 Sharder primal_sharder_;
121 Sharder dual_sharder_;
122};
123
124} // namespace operations_research::pdlp
125
126#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
ShardedQuadraticProgram(QuadraticProgram qp, int num_threads, int num_shards, operations_research::SolverLogger *logger=nullptr)
const Sharder & DualSharder() const
Returns a Sharder intended for dual vectors.
ShardedQuadraticProgram & operator=(const ShardedQuadraticProgram &)=delete
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)
double upper_bound
double lower_bound
int constraint_index
Validation utilities for solvers.proto.