Google OR-Tools v9.11
a fast and portable software suite for combinatorial optimization
Loading...
Searching...
No Matches
sharder.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_SHARDER_H_
15#define PDLP_SHARDER_H_
16
17#include <cstdint>
18#include <functional>
19#include <type_traits>
20#include <vector>
21
22#include "Eigen/Core"
23#include "Eigen/SparseCore"
24#include "absl/log/check.h"
26
28
29// This class represents a way to shard elements for multi-threading. The
30// elements may be entries of a vector or the columns of a (column-major)
31// matrix. The shards are selected to have roughly the same mass, where the mass
32// of an entry depends on the constructor used. See the free functions below and
33// in the .cc file for example usage.
34class Sharder {
35 public:
36 // These are public aliases for convenience. They will change only if there
37 // are breaking changes in Eigen.
38 using ConstSparseColumnBlock = ::Eigen::Block<
39 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>,
40 /*BlockRows=*/Eigen::Dynamic, /*BlockCols=*/Eigen::Dynamic,
41 /*InnerPanel=*/true>;
43 ::Eigen::Block<Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>,
44 /*BlockRows=*/Eigen::Dynamic, /*BlockCols=*/Eigen::Dynamic,
45 /*InnerPanel=*/true>;
46
47 // This class extracts a particular shard of vectors or matrices passed to it.
48 // See `ParallelForEachShard()`.
49 // Caution: Like `absl::Span`, `Shard::operator()` returns mutable or
50 // immutable views into the vector or matrix argument. The underlying object
51 // must outlive the view.
52 // Extra Caution: The const& arguments for the immutable views can bind to
53 // temporary objects, e.g., shard(3*a) will create a view into the "3*a"
54 // object that will be destroyed immediately after the shard is created.
55 class Shard {
56 public:
57 // Returns this shard of `vector`.
58 Eigen::VectorBlock<const Eigen::VectorXd> operator()(
59 const Eigen::VectorXd& vector) const {
60 CHECK_EQ(vector.size(), parent_.NumElements());
61 return vector.segment(parent_.ShardStart(shard_num_),
62 parent_.ShardSize(shard_num_));
63 }
64 // Returns this shard of `vector` in mutable form.
65 Eigen::VectorBlock<Eigen::VectorXd> operator()(
66 Eigen::VectorXd& vector) const {
67 CHECK_EQ(vector.size(), parent_.NumElements());
68 return vector.segment(parent_.ShardStart(shard_num_),
69 parent_.ShardSize(shard_num_));
70 }
71 // Returns this shard of `vector`.
72 Eigen::VectorBlock<const Eigen::VectorXd> operator()(
73 Eigen::VectorBlock<const Eigen::VectorXd> vector) const {
74 CHECK_EQ(vector.size(), parent_.NumElements());
75 return Eigen::VectorBlock<const Eigen::VectorXd>(
76 vector.nestedExpression(),
77 vector.startRow() + parent_.ShardStart(shard_num_),
78 parent_.ShardSize(shard_num_));
79 }
80 // Returns this shard of `vector` in mutable form.
81 Eigen::VectorBlock<Eigen::VectorXd> operator()(
82 Eigen::VectorBlock<Eigen::VectorXd> vector) const {
83 CHECK_EQ(vector.size(), parent_.NumElements());
84 return Eigen::VectorBlock<Eigen::VectorXd>(
85 vector.nestedExpression(),
86 vector.startRow() + parent_.ShardStart(shard_num_),
87 parent_.ShardSize(shard_num_));
88 }
89 // Returns this shard of `diag`. Note that the shard is a *square* diagonal
90 // matrix, not a block of columns of original length.
91 auto operator()(const Eigen::DiagonalMatrix<double, Eigen::Dynamic>& diag)
92 const -> decltype(diag.diagonal().segment(0, 0).asDiagonal()) {
93 CHECK_EQ(diag.diagonal().size(), parent_.NumElements());
94 return diag.diagonal()
95 .segment(parent_.ShardStart(shard_num_),
96 parent_.ShardSize(shard_num_))
97 .asDiagonal();
98 }
99 // Returns this shard of the columns of `matrix`.
101 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix)
102 const {
103 CHECK_EQ(matrix.cols(), parent_.NumElements());
104 auto result = matrix.middleCols(parent_.ShardStart(shard_num_),
105 parent_.ShardSize(shard_num_));
106 // This is a guard against implicit conversions, because the return type
107 // of `middleCols` is not 100% clear from the documentation.
108 static_assert(
109 std::is_same<decltype(result), ConstSparseColumnBlock>::value,
110 "The return type of middleCols changed!");
111 return result;
112 }
113 // Returns this shard of the columns of `matrix` in mutable form.
115 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix) const {
116 CHECK_EQ(matrix.cols(), parent_.NumElements());
117 auto result = matrix.middleCols(parent_.ShardStart(shard_num_),
118 parent_.ShardSize(shard_num_));
119 // This is a guard against implicit conversions, because the return type
120 // of `middleCols` is not 100% clear from the documentation.
121 static_assert(std::is_same<decltype(result), SparseColumnBlock>::value,
122 "The return type of middleCols changed!");
123 return result;
124 }
125 // A non-negative identifier for this shard, less than `NumShards()` of the
126 // parent `Sharder`.
127 int Index() const { return shard_num_; }
128
129 private:
130 friend class Sharder;
131 Shard(int shard_num, const Sharder* parent)
132 : shard_num_(shard_num), parent_(*parent) {
133 CHECK_NE(parent, nullptr);
134 CHECK_GE(shard_num, 0);
135 CHECK_LT(shard_num, parent->NumShards());
136 }
137 const int shard_num_;
138 const Sharder& parent_;
139 };
140
141 // Creates a `Sharder` for problems with `num_elements` elements and mass of
142 // each element given by `element_mass`. Each shard will have roughly the same
143 // mass. The number of shards in the resulting `Sharder` will be approximately
144 // `num_shards` but may differ. The `thread_pool` will be used for parallel
145 // operations executed by e.g. `ParallelForEachShard()`. The `thread_pool` may
146 // be nullptr, which means work will be executed in the same thread. If
147 // `thread_pool` is not nullptr, the underlying object is not owned and must
148 // outlive the `Sharder`.
149 Sharder(int64_t num_elements, int num_shards, ThreadPool* thread_pool,
150 const std::function<int64_t(int64_t)>& element_mass);
151
152 // Creates a `Sharder` for problems with `num_elements` elements and unit
153 // mass. This constructor exploits having all element mass equal to 1 to take
154 // time proportional to `num_shards` instead of `num_elements`. Also see the
155 // comments above the first constructor.
156 Sharder(int64_t num_elements, int num_shards, ThreadPool* thread_pool);
157
158 // Creates a `Sharder` for processing `matrix`. The elements correspond to
159 // columns of `matrix` and have mass linear in the number of non-zeros. Also
160 // see the comments above the first constructor.
161 Sharder(const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
162 int num_shards, ThreadPool* thread_pool)
163 : Sharder(matrix.cols(), num_shards, thread_pool, [&matrix](int64_t col) {
164 return 1 + 1 * matrix.col(col).nonZeros();
165 }) {}
166
167 // Constructs a `Sharder` with the same thread pool as `other_sharder`, for
168 // problems with `num_elements` elements and unit mass. The number of shards
169 // will be approximately the same as that of `other_sharder`. Also see the
170 // comments on the first constructor.
171 Sharder(const Sharder& other_sharder, int64_t num_elements);
172
173 // `Sharder` may be moved, but not copied.
174 // Moved-from objects may be in an invalid state. The only methods that may be
175 // called on a moved-from object are the destructor or `operator=`.
176 Sharder(const Sharder& other) = delete;
177 Sharder& operator=(const Sharder& other) = delete;
178 Sharder(Sharder&& other) = default;
179 Sharder& operator=(Sharder&& other) = default;
180
181 int NumShards() const { return static_cast<int>(shard_starts_.size()) - 1; }
182
183 // The number of elements that are split into shards.
184 int64_t NumElements() const { return shard_starts_.back(); }
185
186 int64_t ShardSize(int shard) const {
187 CHECK_GE(shard, 0);
188 CHECK_LT(shard, NumShards());
189 return shard_starts_[shard + 1] - shard_starts_[shard];
190 }
191
192 int64_t ShardStart(int shard) const {
193 CHECK_GE(shard, 0);
194 CHECK_LT(shard, NumShards());
195 return shard_starts_[shard];
196 }
197
198 int64_t ShardMass(int shard) const {
199 CHECK_GE(shard, 0);
200 CHECK_LT(shard, NumShards());
201 return shard_masses_[shard];
202 }
203
204 // Runs `func` on each of the shards.
206 const std::function<void(const Shard&)>& func) const;
207
208 // Runs `func` on each of the shards and sums the results.
210 const std::function<double(const Shard&)>& func) const;
211
212 // Runs `func` on each of the shards. Returns true iff all shards returned
213 // true.
215 const std::function<bool(const Shard&)>& func) const;
216
217 // Public for testing only.
218 const std::vector<int64_t>& ShardStartsForTesting() const {
219 return shard_starts_;
220 }
221
222 private:
223 // Size: `NumShards() + 1`. The first entry is 0 and the last entry is
224 // `NumElements()`. The entries are sorted in increasing order and are unique.
225 // Note that {0} is valid and indicates zero elements split into zero shards.
226 std::vector<int64_t> shard_starts_;
227 // Size: `NumShards()`. The mass of each shard.
228 std::vector<int64_t> shard_masses_;
229 // NOT owned. May be nullptr.
230 ThreadPool* thread_pool_;
231};
232
233// Like `matrix.transpose() * vector` but executed in parallel using `sharder`.
234// The size of `sharder` must match the number of columns in `matrix`. To ensure
235// good parallelization `matrix` should have (roughly) the same location of
236// non-zeros as the `matrix` used when constructing `sharder`.
238 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
239 const Eigen::VectorXd& vector, const Sharder& sharder);
240
242// The following functions use `sharder` to compute a vector operation in
243// parallel. `sharder` should have the same size as the vector(s). For best
244// performance `sharder` should have been created with the `Sharder(int64_t,
245// int, ThreadPool*)` constructor.
247
248// Like `dest.setZero(sharder.NumElements())`. Note that if `dest.size() !=
249// sharder.NumElements()`, `dest` will be resized.
250void SetZero(const Sharder& sharder, Eigen::VectorXd& dest);
251
252// Like `VectorXd::Zero(sharder.NumElements())`.
253Eigen::VectorXd ZeroVector(const Sharder& sharder);
254
255// Like `VectorXd::Ones(sharder.NumElements())`.
256Eigen::VectorXd OnesVector(const Sharder& sharder);
257
258// Like `dest += scale * increment`.
259void AddScaledVector(double scale, const Eigen::VectorXd& increment,
260 const Sharder& sharder, Eigen::VectorXd& dest);
262// Like `dest = vec`. `dest` is resized if needed.
263void AssignVector(const Eigen::VectorXd& vec, const Sharder& sharder,
264 Eigen::VectorXd& dest);
266// Returns a copy of `vec`.
267Eigen::VectorXd CloneVector(const Eigen::VectorXd& vec, const Sharder& sharder);
268
269// Like `dest = dest.cwiseProduct(scale)`.
270void CoefficientWiseProductInPlace(const Eigen::VectorXd& scale,
271 const Sharder& sharder,
272 Eigen::VectorXd& dest);
273
274// Like `dest = dest.cwiseQuotient(scale)`.
275void CoefficientWiseQuotientInPlace(const Eigen::VectorXd& scale,
276 const Sharder& sharder,
277 Eigen::VectorXd& dest);
278
279// Like `v1.dot(v2)`.
280double Dot(const Eigen::VectorXd& v1, const Eigen::VectorXd& v2,
281 const Sharder& sharder);
282// Like `vector.lpNorm<Eigen::Infinity>()`, a.k.a. LInf norm.
283double LInfNorm(const Eigen::VectorXd& vector, const Sharder& sharder);
284// Like `vector.lpNorm<1>()`, a.k.a. L_1 norm.
285double L1Norm(const Eigen::VectorXd& vector, const Sharder& sharder);
286// Like `vector.squaredNorm()`.
287double SquaredNorm(const Eigen::VectorXd& vector, const Sharder& sharder);
288// Like `vector.norm()`.
289double Norm(const Eigen::VectorXd& vector, const Sharder& sharder);
290
291// Like `(vector1 - vector2).squaredNorm()`.
292double SquaredDistance(const Eigen::VectorXd& vector1,
293 const Eigen::VectorXd& vector2, const Sharder& sharder);
294// Like `(vector1 - vector2).norm()`.
295double Distance(const Eigen::VectorXd& vector1, const Eigen::VectorXd& vector2,
296 const Sharder& sharder);
298// `ScaledL1Norm` is omitted because it's not needed (yet).
299
300// LInf norm of a rescaled vector, i.e.,
301// `vector.cwiseProduct(scale).lpNorm<Eigen::Infinity>()`.
302double ScaledLInfNorm(const Eigen::VectorXd& vector,
303 const Eigen::VectorXd& scale, const Sharder& sharder);
304// Squared L2 norm of a rescaled vector, i.e.,
305// `vector.cwiseProduct(scale).squaredNorm()`.
306double ScaledSquaredNorm(const Eigen::VectorXd& vector,
307 const Eigen::VectorXd& scale, const Sharder& sharder);
308// L2 norm of a rescaled vector, i.e., `vector.cwiseProduct(scale).norm()`.
309double ScaledNorm(const Eigen::VectorXd& vector, const Eigen::VectorXd& scale,
310 const Sharder& sharder);
313// The functions below compute norms of the columns of a scaled matrix. The
314// (i,j) entry of the scaled matrix equals `matrix[i,j] * row_scaling_vec[i]
315// * col_scaling_vec[j]`. To ensure good parallelization `matrix` should have
316// (roughly) the same location of non-zeros as the `matrix` used to construct
317// `sharder`. The size of `sharder` must match the number of columns in
318// `matrix`.
320
321// Computes the LInf norm of each column of a scaled `matrix`.
322Eigen::VectorXd ScaledColLInfNorm(
323 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
324 const Eigen::VectorXd& row_scaling_vec,
325 const Eigen::VectorXd& col_scaling_vec, const Sharder& sharder);
326// Computes the L2 norm of each column of a scaled `matrix`.
327Eigen::VectorXd ScaledColL2Norm(
328 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
329 const Eigen::VectorXd& row_scaling_vec,
330 const Eigen::VectorXd& col_scaling_vec, const Sharder& sharder);
332} // namespace operations_research::pdlp
333
334#endif // PDLP_SHARDER_H_
Eigen::VectorBlock< const Eigen::VectorXd > operator()(Eigen::VectorBlock< const Eigen::VectorXd > vector) const
Returns this shard of vector.
Definition sharder.h:72
SparseColumnBlock operator()(Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix) const
Returns this shard of the columns of matrix in mutable form.
Definition sharder.h:114
ConstSparseColumnBlock operator()(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix) const
Returns this shard of the columns of matrix.
Definition sharder.h:100
Eigen::VectorBlock< Eigen::VectorXd > operator()(Eigen::VectorXd &vector) const
Returns this shard of vector in mutable form.
Definition sharder.h:65
Eigen::VectorBlock< const Eigen::VectorXd > operator()(const Eigen::VectorXd &vector) const
Returns this shard of vector.
Definition sharder.h:58
Eigen::VectorBlock< Eigen::VectorXd > operator()(Eigen::VectorBlock< Eigen::VectorXd > vector) const
Returns this shard of vector in mutable form.
Definition sharder.h:81
auto operator()(const Eigen::DiagonalMatrix< double, Eigen::Dynamic > &diag) const -> decltype(diag.diagonal().segment(0, 0).asDiagonal())
Definition sharder.h:91
const std::vector< int64_t > & ShardStartsForTesting() const
Public for testing only.
Definition sharder.h:218
Sharder(const Sharder &other)=delete
Sharder(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, int num_shards, ThreadPool *thread_pool)
Definition sharder.h:161
Sharder & operator=(Sharder &&other)=default
bool ParallelTrueForAllShards(const std::function< bool(const Shard &)> &func) const
Definition sharder.cc:152
double ParallelSumOverShards(const std::function< double(const Shard &)> &func) const
Runs func on each of the shards and sums the results.
Definition sharder.cc:143
Sharder & operator=(const Sharder &other)=delete
int64_t ShardSize(int shard) const
Definition sharder.h:186
::Eigen::Block< Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t >, Eigen::Dynamic, Eigen::Dynamic, true > SparseColumnBlock
Definition sharder.h:42
::Eigen::Block< const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t >, Eigen::Dynamic, Eigen::Dynamic, true > ConstSparseColumnBlock
Definition sharder.h:38
Sharder(Sharder &&other)=default
int64_t ShardMass(int shard) const
Definition sharder.h:198
Sharder(int64_t num_elements, int num_shards, ThreadPool *thread_pool, const std::function< int64_t(int64_t)> &element_mass)
Definition sharder.cc:36
void ParallelForEachShard(const std::function< void(const Shard &)> &func) const
Runs func on each of the shards.
Definition sharder.cc:109
int64_t NumElements() const
The number of elements that are split into shards.
Definition sharder.h:184
int64_t ShardStart(int shard) const
Definition sharder.h:192
int64_t value
ColIndex col
Definition markowitz.cc:187
Validation utilities for solvers.proto.
double SquaredNorm(const VectorXd &vector, const Sharder &sharder)
Definition sharder.cc:248
void SetZero(const Sharder &sharder, VectorXd &dest)
Definition sharder.cc:178
double ScaledNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition sharder.cc:286
double Dot(const VectorXd &v1, const VectorXd &v2, const Sharder &sharder)
Definition sharder.cc:230
double SquaredDistance(const VectorXd &vector1, const VectorXd &vector2, const Sharder &sharder)
Definition sharder.cc:257
double LInfNorm(const VectorXd &vector, const Sharder &sharder)
Definition sharder.cc:235
double Distance(const VectorXd &vector1, const VectorXd &vector2, const Sharder &sharder)
Definition sharder.cc:264
VectorXd TransposedMatrixVectorProduct(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &vector, const Sharder &sharder)
Definition sharder.cc:163
double ScaledLInfNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition sharder.cc:269
double ScaledSquaredNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition sharder.cc:279
VectorXd ScaledColLInfNorm(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &row_scaling_vec, const VectorXd &col_scaling_vec, const Sharder &sharder)
Definition sharder.cc:291
void AddScaledVector(const double scale, const VectorXd &increment, const Sharder &sharder, VectorXd &dest)
Definition sharder.cc:197
void CoefficientWiseProductInPlace(const VectorXd &scale, const Sharder &sharder, VectorXd &dest)
Definition sharder.cc:216
void CoefficientWiseQuotientInPlace(const VectorXd &scale, const Sharder &sharder, VectorXd &dest)
Definition sharder.cc:223
VectorXd ZeroVector(const Sharder &sharder)
Like VectorXd::Zero(sharder.NumElements()).
Definition sharder.cc:184
VectorXd ScaledColL2Norm(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &row_scaling_vec, const VectorXd &col_scaling_vec, const Sharder &sharder)
Definition sharder.cc:313
double L1Norm(const VectorXd &vector, const Sharder &sharder)
Definition sharder.cc:243
VectorXd CloneVector(const VectorXd &vec, const Sharder &sharder)
Definition sharder.cc:210
VectorXd OnesVector(const Sharder &sharder)
Like VectorXd::Ones(sharder.NumElements()).
Definition sharder.cc:190
double Norm(const VectorXd &vector, const Sharder &sharder)
Definition sharder.cc:253
void AssignVector(const VectorXd &vec, const Sharder &sharder, VectorXd &dest)
Definition sharder.cc:204