OR-Tools  9.3
sharder.cc
Go to the documentation of this file.
1// Copyright 2010-2021 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
15
16#include <algorithm>
17#include <cmath>
18#include <cstdint>
19#include <functional>
20#include <vector>
21
22#include "Eigen/Core"
23#include "Eigen/SparseCore"
24#include "absl/synchronization/blocking_counter.h"
25#include "absl/time/time.h"
29#include "ortools/base/timer.h"
30
32
33using ::Eigen::VectorXd;
34
35Sharder::Sharder(const int64_t num_elements, const int num_shards,
36 ThreadPool* const thread_pool,
37 const std::function<int64_t(int64_t)>& element_mass)
38 : thread_pool_(thread_pool) {
39 CHECK_GE(num_elements, 0);
40 if (num_elements == 0) {
41 shard_starts_.push_back(0);
42 return;
43 }
44 CHECK_GE(num_shards, 1);
45 int64_t overall_mass = 0;
46 for (int64_t elem = 0; elem < num_elements; ++elem) {
47 overall_mass += element_mass(elem);
48 }
49 shard_starts_.push_back(0);
50 int64_t this_shard_mass = element_mass(0);
51 for (int64_t elem = 1; elem < num_elements; ++elem) {
52 int64_t this_elem_mass = element_mass(elem);
53 if (this_shard_mass + (this_elem_mass / 2) >= overall_mass / num_shards) {
54 // this elem starts a new shard
55 shard_masses_.push_back(this_shard_mass);
56 shard_starts_.push_back(elem);
57 this_shard_mass = this_elem_mass;
58 } else {
59 this_shard_mass += this_elem_mass;
60 }
61 }
62 shard_starts_.push_back(num_elements);
63 shard_masses_.push_back(this_shard_mass);
64 CHECK_EQ(NumShards(), shard_masses_.size());
65}
66
67Sharder::Sharder(const int64_t num_elements, const int num_shards,
68 ThreadPool* const thread_pool)
69 : thread_pool_(thread_pool) {
70 CHECK_GE(num_elements, 0);
71 if (num_elements == 0) {
72 shard_starts_.push_back(0);
73 return;
74 }
75 CHECK_GE(num_shards, 1);
76 shard_starts_.reserve(num_shards + 1);
77 shard_masses_.reserve(num_shards);
78 for (int shard = 0; shard < num_shards; ++shard) {
79 const int64_t this_shard_start = ((num_elements * shard) / num_shards);
80 const int64_t next_shard_start =
81 ((num_elements * (shard + 1)) / num_shards);
82 if (next_shard_start - this_shard_start > 0) {
83 shard_starts_.push_back(this_shard_start);
84 shard_masses_.push_back(next_shard_start - this_shard_start);
85 }
86 }
87 shard_starts_.push_back(num_elements);
88 CHECK_EQ(NumShards(), shard_masses_.size());
89}
90
91Sharder::Sharder(const Sharder& other_sharder, const int64_t num_elements)
92 // The std::max() protects against other_sharder.NumShards() == 0, which
93 // will happen if other_sharder had num_elements == 0.
94 : Sharder(num_elements, std::max(1, other_sharder.NumShards()),
95 other_sharder.thread_pool_) {}
96
98 const std::function<void(const Shard&)>& func) const {
99 if (thread_pool_) {
100 absl::BlockingCounter counter(NumShards());
101 VLOG(2) << "Starting ParallelForEachShard()";
102 for (int shard_num = 0; shard_num < NumShards(); ++shard_num) {
103 thread_pool_->Schedule([&, shard_num]() {
104 WallTimer timer;
105 if (VLOG_IS_ON(2)) {
106 timer.Start();
107 }
108 func(Shard(shard_num, this));
109 if (VLOG_IS_ON(2)) {
110 timer.Stop();
111 VLOG(2) << "Shard " << shard_num << " with " << ShardSize(shard_num)
112 << " elements and " << ShardMass(shard_num)
113 << " mass finished with "
114 << ShardMass(shard_num) /
115 std::max(int64_t{1}, absl::ToInt64Microseconds(
116 timer.GetDuration()))
117 << " mass/usec.";
118 }
119 counter.DecrementCount();
120 });
121 }
122 counter.Wait();
123 VLOG(2) << "Done ParallelForEachShard()";
124 } else {
125 for (int shard_num = 0; shard_num < NumShards(); ++shard_num) {
126 func(Shard(shard_num, this));
127 }
128 }
129}
130
132 const std::function<double(const Shard&)>& func) const {
133 VectorXd local_sums(NumShards());
134 ParallelForEachShard([&](const Sharder::Shard& shard) {
135 local_sums[shard.Index()] = func(shard);
136 });
137 return local_sums.sum();
138}
139
141 const std::function<bool(const Shard&)>& func) const {
142 // Recall std::vector<bool> is not thread-safe.
143 std::vector<int> local_result(NumShards());
144 ParallelForEachShard([&](const Sharder::Shard& shard) {
145 local_result[shard.Index()] = static_cast<int>(func(shard));
146 });
147 return std::all_of(local_result.begin(), local_result.end(),
148 [](const int v) { return static_cast<bool>(v); });
149}
150
152 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
153 const VectorXd& vector, const Sharder& sharder) {
154 CHECK_EQ(vector.size(), matrix.rows());
155 VectorXd answer(matrix.cols());
156 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
157 shard(answer) = shard(matrix).transpose() * vector;
158 });
159 return answer;
160}
161
162void AddScaledVector(const double scale, const VectorXd& increment,
163 const Sharder& sharder, VectorXd& dest) {
164 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
165 shard(dest) += scale * shard(increment);
166 });
167}
168
169void AssignVector(const VectorXd& vec, const Sharder& sharder, VectorXd& dest) {
170 dest.resize(vec.size());
171 sharder.ParallelForEachShard(
172 [&](const Sharder::Shard& shard) { shard(dest) = shard(vec); });
173}
174
175VectorXd CloneVector(const VectorXd& vec, const Sharder& sharder) {
176 VectorXd dest;
177 AssignVector(vec, sharder, dest);
178 return dest;
179}
180
181// Like vector = vector.cwiseProduct(scale).
182void CoefficientWiseProductInPlace(const VectorXd& scale,
183 const Sharder& sharder, VectorXd& dest) {
184 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
185 shard(dest) = shard(dest).cwiseProduct(shard(scale));
186 });
187}
188
189// Like vector = vector.cwiseQuotient(scale).
190void CoefficientWiseQuotientInPlace(const VectorXd& scale,
191 const Sharder& sharder, VectorXd& dest) {
192 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
193 shard(dest) = shard(dest).cwiseQuotient(shard(scale));
194 });
195}
196
197double Dot(const VectorXd& v1, const VectorXd& v2, const Sharder& sharder) {
198 return sharder.ParallelSumOverShards(
199 [&](const Sharder::Shard& shard) { return shard(v1).dot(shard(v2)); });
200}
201
202double LInfNorm(const VectorXd& vector, const Sharder& sharder) {
203 VectorXd local_max(sharder.NumShards());
204 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
205 local_max[shard.Index()] = shard(vector).lpNorm<Eigen::Infinity>();
206 });
207 return local_max.lpNorm<Eigen::Infinity>();
208}
209
210double L1Norm(const VectorXd& vector, const Sharder& sharder) {
211 return sharder.ParallelSumOverShards(
212 [&](const Sharder::Shard& shard) { return shard(vector).lpNorm<1>(); });
213}
214
215double SquaredNorm(const VectorXd& vector, const Sharder& sharder) {
216 return sharder.ParallelSumOverShards(
217 [&](const Sharder::Shard& shard) { return shard(vector).squaredNorm(); });
218}
219
220double Norm(const VectorXd& vector, const Sharder& sharder) {
221 return std::sqrt(SquaredNorm(vector, sharder));
222}
223
224double SquaredDistance(const VectorXd& vector1, const VectorXd& vector2,
225 const Sharder& sharder) {
226 return sharder.ParallelSumOverShards([&](const Sharder::Shard& shard) {
227 return (shard(vector1) - shard(vector2)).squaredNorm();
228 });
229}
230
231double Distance(const VectorXd& vector1, const VectorXd& vector2,
232 const Sharder& sharder) {
233 return std::sqrt(SquaredDistance(vector1, vector2, sharder));
234}
235
236double ScaledLInfNorm(const VectorXd& vector, const VectorXd& scale,
237 const Sharder& sharder) {
238 VectorXd local_max(sharder.NumShards());
239 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
240 local_max[shard.Index()] =
241 shard(vector).cwiseProduct(shard(scale)).lpNorm<Eigen::Infinity>();
242 });
243 return local_max.lpNorm<Eigen::Infinity>();
244}
245
246double ScaledSquaredNorm(const VectorXd& vector, const VectorXd& scale,
247 const Sharder& sharder) {
248 return sharder.ParallelSumOverShards([&](const Sharder::Shard& shard) {
249 return shard(vector).cwiseProduct(shard(scale)).squaredNorm();
250 });
251}
252
253double ScaledNorm(const VectorXd& vector, const VectorXd& scale,
254 const Sharder& sharder) {
255 return std::sqrt(ScaledSquaredNorm(vector, scale, sharder));
256}
257
259 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
260 const VectorXd& row_scaling_vec, const VectorXd& col_scaling_vec,
261 const Sharder& sharder) {
262 CHECK_EQ(matrix.cols(), col_scaling_vec.size());
263 CHECK_EQ(matrix.rows(), row_scaling_vec.size());
264 VectorXd answer(matrix.cols());
265 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
266 auto matrix_shard = shard(matrix);
267 auto col_scaling_shard = shard(col_scaling_vec);
268 for (int64_t col_num = 0; col_num < shard(matrix).outerSize(); ++col_num) {
269 double max = 0.0;
270 for (decltype(matrix_shard)::InnerIterator it(matrix_shard, col_num); it;
271 ++it) {
272 max = std::max(max, std::abs(it.value() * row_scaling_vec[it.row()]));
273 }
274 shard(answer)[col_num] = max * std::abs(col_scaling_shard[col_num]);
275 }
276 });
277 return answer;
278}
279
281 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
282 const VectorXd& row_scaling_vec, const VectorXd& col_scaling_vec,
283 const Sharder& sharder) {
284 CHECK_EQ(matrix.cols(), col_scaling_vec.size());
285 CHECK_EQ(matrix.rows(), row_scaling_vec.size());
286 VectorXd answer(matrix.cols());
287 sharder.ParallelForEachShard([&](const Sharder::Shard& shard) {
288 auto matrix_shard = shard(matrix);
289 auto col_scaling_shard = shard(col_scaling_vec);
290 for (int64_t col_num = 0; col_num < shard(matrix).outerSize(); ++col_num) {
291 double sum_of_squares = 0.0;
292 for (decltype(matrix_shard)::InnerIterator it(matrix_shard, col_num); it;
293 ++it) {
294 sum_of_squares +=
295 MathUtil::Square(it.value() * row_scaling_vec[it.row()]);
296 }
297 shard(answer)[col_num] =
298 std::sqrt(sum_of_squares) * std::abs(col_scaling_shard[col_num]);
299 }
300 });
301 return answer;
302}
303
305 const Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t>& matrix,
306 const Sharder& sharder) {
307 return sharder.ParallelTrueForAllShards([&](const Sharder::Shard& shard) {
308 auto matrix_shard = shard(matrix);
309 const int64_t col_offset = sharder.ShardStart(shard.Index());
310 for (int64_t col_idx = 0; col_idx < matrix_shard.outerSize(); ++col_idx) {
311 for (decltype(matrix_shard)::InnerIterator it(matrix_shard, col_idx); it;
312 ++it) {
313 if (it.row() != (col_offset + it.col())) return false;
314 }
315 }
316 return true;
317 });
318}
319
320} // namespace operations_research::pdlp
int64_t max
Definition: alldiff_cst.cc:140
#define CHECK_EQ(val1, val2)
Definition: base/logging.h:703
#define CHECK_GE(val1, val2)
Definition: base/logging.h:707
#define VLOG(verboselevel)
Definition: base/logging.h:984
void Start()
Definition: timer.h:31
void Stop()
Definition: timer.h:39
absl::Duration GetDuration() const
Definition: timer.h:48
static T Square(const T x)
Definition: mathutil.h:101
void Schedule(std::function< void()> closure)
Definition: threadpool.cc:77
Sharder(int64_t num_elements, int num_shards, ThreadPool *thread_pool, const std::function< int64_t(int64_t)> &element_mass)
Definition: sharder.cc:35
double ParallelSumOverShards(const std::function< double(const Shard &)> &func) const
Definition: sharder.cc:131
void ParallelForEachShard(const std::function< void(const Shard &)> &func) const
Definition: sharder.cc:97
bool ParallelTrueForAllShards(const std::function< bool(const Shard &)> &func) const
Definition: sharder.cc:140
int64_t ShardSize(int shard) const
Definition: sharder.h:186
int64_t ShardStart(int shard) const
Definition: sharder.h:192
int64_t ShardMass(int shard) const
Definition: sharder.h:198
double SquaredNorm(const VectorXd &vector, const Sharder &sharder)
Definition: sharder.cc:215
double ScaledNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition: sharder.cc:253
double Dot(const VectorXd &v1, const VectorXd &v2, const Sharder &sharder)
Definition: sharder.cc:197
double SquaredDistance(const VectorXd &vector1, const VectorXd &vector2, const Sharder &sharder)
Definition: sharder.cc:224
double LInfNorm(const VectorXd &vector, const Sharder &sharder)
Definition: sharder.cc:202
double Distance(const VectorXd &vector1, const VectorXd &vector2, const Sharder &sharder)
Definition: sharder.cc:231
VectorXd TransposedMatrixVectorProduct(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &vector, const Sharder &sharder)
Definition: sharder.cc:151
double ScaledLInfNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition: sharder.cc:236
double ScaledSquaredNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
Definition: sharder.cc:246
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:258
void AddScaledVector(const double scale, const VectorXd &increment, const Sharder &sharder, VectorXd &dest)
Definition: sharder.cc:162
void CoefficientWiseProductInPlace(const VectorXd &scale, const Sharder &sharder, VectorXd &dest)
Definition: sharder.cc:182
void CoefficientWiseQuotientInPlace(const VectorXd &scale, const Sharder &sharder, VectorXd &dest)
Definition: sharder.cc:190
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:280
double L1Norm(const VectorXd &vector, const Sharder &sharder)
Definition: sharder.cc:210
VectorXd CloneVector(const VectorXd &vec, const Sharder &sharder)
Definition: sharder.cc:175
bool IsDiagonal(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const Sharder &sharder)
Definition: sharder.cc:304
double Norm(const VectorXd &vector, const Sharder &sharder)
Definition: sharder.cc:220
void AssignVector(const VectorXd &vec, const Sharder &sharder, VectorXd &dest)
Definition: sharder.cc:169
STL namespace.
#define VLOG_IS_ON(verboselevel)
Definition: vlog_is_on.h:44