24#include "Eigen/SparseCore"
25#include "absl/random/distributions.h"
26#include "gmock/gmock.h"
27#include "gtest/gtest.h"
35using ::Eigen::DiagonalMatrix;
36using ::Eigen::VectorXd;
37using ::testing::DoubleNear;
38using ::testing::ElementsAre;
40using Shard = Sharder::Shard;
46Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> TestSparseMatrix() {
47 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat(3, 4);
48 mat.coeffRef(0, 0) = 7;
49 mat.coeffRef(0, 1) = -0.5;
50 mat.coeffRef(1, 0) = 1;
51 mat.coeffRef(1, 2) = 3;
52 mat.coeffRef(1, 3) = 2;
53 mat.coeffRef(2, 0) = -1;
54 mat.coeffRef(2, 3) = 5;
61Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> LargeSparseMatrix(
64 std::mt19937 rand(48709241);
65 std::vector<Eigen::Triplet<double, int64_t>> triplets;
66 for (int64_t
col = 0;
col < size; ++
col) {
69 row += absl::Uniform(rand, 1,
col + 2);
71 double value = absl::Uniform(rand, 1, 10);
76 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat(size, size);
77 mat.setFromTriplets(triplets.begin(), triplets.end());
83void VerifySharder(
const Sharder& sharder,
const int target_num_shards,
84 const std::vector<int64_t>& element_masses) {
85 int64_t num_elements = element_masses.size();
86 int num_shards = sharder.NumShards();
87 ASSERT_EQ(sharder.NumElements(), num_elements);
88 ASSERT_GE(num_elements, 1);
89 ASSERT_GE(num_shards, 1);
90 int64_t elements_so_far = 0;
91 for (
int shard = 0; shard < num_shards; ++shard) {
92 int64_t shard_start = sharder.ShardStart(shard);
93 EXPECT_EQ(shard_start, elements_so_far) <<
" in shard: " << shard;
94 int64_t shard_mass = 0;
95 EXPECT_GE(sharder.ShardSize(shard), 1) <<
" in shard: " << shard;
96 EXPECT_GE(sharder.ShardMass(shard), 1) <<
" in shard: " << shard;
97 for (int64_t i = 0; i < sharder.ShardSize(shard); ++i) {
98 shard_mass += element_masses[shard_start + i];
100 EXPECT_EQ(shard_mass, sharder.ShardMass(shard)) <<
" in shard: " << shard;
101 elements_so_far += sharder.ShardSize(shard);
103 EXPECT_EQ(elements_so_far, num_elements);
105 EXPECT_LE(num_shards, 2 * target_num_shards);
106 ASSERT_GE(target_num_shards, 1);
107 const int64_t overall_mass =
108 std::accumulate(element_masses.begin(), element_masses.end(), int64_t{0});
109 const int64_t max_element_mass =
110 *std::max_element(element_masses.begin(), element_masses.end());
111 const int64_t upper_mass_limit =
std::max(
115 const int64_t lower_mass_limit =
116 overall_mass / target_num_shards -
118 for (
int shard = 0; shard < sharder.NumShards(); ++shard) {
119 EXPECT_LE(sharder.ShardMass(shard), upper_mass_limit)
120 <<
" in shard: " << shard;
121 if (shard + 1 < sharder.NumShards()) {
122 EXPECT_GE(sharder.ShardMass(shard), lower_mass_limit)
123 <<
" in shard: " << shard;
128TEST(SharderTest, SharderFromMatrix) {
129 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat =
131 Sharder sharder(mat, 2,
nullptr);
132 VerifySharder(sharder, 2, {4, 2, 2, 3});
135TEST(SharderTest, UniformSharder) {
136 Sharder sharder(10, 3,
nullptr);
137 VerifySharder(sharder, 3, {1, 1, 1, 1, 1, 1, 1, 1, 1, 1});
140TEST(SharderTest, UniformSharderFromOtherSharder) {
141 Sharder other_sharder(5, 3,
nullptr);
142 Sharder sharder(other_sharder, 10);
143 VerifySharder(sharder, other_sharder.NumShards(),
144 {1, 1, 1, 1, 1, 1, 1, 1, 1, 1});
147TEST(SharderTest, UniformSharderExcessiveShards) {
148 Sharder sharder(5, 7,
nullptr);
149 EXPECT_THAT(sharder.ShardStartsForTesting(), ElementsAre(0, 1, 2, 3, 4, 5));
150 VerifySharder(sharder, 7, {1, 1, 1, 1, 1});
153TEST(SharderTest, UniformSharderOneShard) {
154 Sharder sharder(5, 1,
nullptr);
155 EXPECT_THAT(sharder.ShardStartsForTesting(), ElementsAre(0, 5));
156 VerifySharder(sharder, 1, {1, 1, 1, 1, 1});
159TEST(SharderTest, UniformSharderOneElementVector) {
160 Sharder sharder(1, 5,
nullptr);
161 EXPECT_THAT(sharder.ShardStartsForTesting(), ElementsAre(0, 1));
162 VerifySharder(sharder, 5, {1});
165TEST(SharderTest, UniformSharderZeroElementVector) {
166 Sharder sharder(0, 3,
nullptr);
167 EXPECT_THAT(sharder.ShardStartsForTesting(), ElementsAre(0));
168 EXPECT_EQ(sharder.NumShards(), 0);
169 EXPECT_EQ(sharder.NumElements(), 0);
170 sharder.ParallelForEachShard([](
const Shard& ) {
171 LOG(
FATAL) <<
"There are no shards so this shouldn't be called.";
175TEST(SharderTest, UniformSharderFromOtherZeroElementSharder) {
176 Sharder empty_sharder(0, 3,
nullptr);
177 EXPECT_THAT(empty_sharder.ShardStartsForTesting(), ElementsAre(0));
178 EXPECT_EQ(empty_sharder.NumShards(), 0);
179 EXPECT_EQ(empty_sharder.NumElements(), 0);
180 Sharder sharder(empty_sharder, 5);
181 EXPECT_THAT(sharder.ShardStartsForTesting(), ElementsAre(0, 5));
182 VerifySharder(sharder, 1, {1, 1, 1, 1, 1});
185TEST(ParallelSumOverShards, SmallExample) {
188 Sharder sharder(vec.size(), 2,
nullptr);
189 const double sum = sharder.ParallelSumOverShards(
190 [&vec](
const Shard& shard) {
return shard(vec).sum(); });
194TEST(ParallelSumOverShards, SmallExampleUsingVectorBlock) {
197 auto vec_block = vec.segment(1, 2);
198 Sharder sharder(vec_block.size(), 2,
nullptr);
199 const double sum = sharder.ParallelSumOverShards(
200 [&vec_block](
const Shard& shard) {
return shard(vec_block).sum(); });
204TEST(ParallelSumOverShards, SmallExampleUsingConstVectorBlock) {
207 const VectorXd& const_vec = vec;
208 auto vec_block = const_vec.segment(1, 2);
209 Sharder sharder(vec_block.size(), 2,
nullptr);
210 const double sum = sharder.ParallelSumOverShards(
211 [&vec_block](
const Shard& shard) {
return shard(vec_block).sum(); });
215TEST(ParallelSumOverShards, SmallExampleUsingDiagonalMatrix) {
216 DiagonalMatrix<double, Eigen::Dynamic> diag{{1, 2, 3}};
217 Sharder sharder(diag.cols(), 2,
nullptr);
218 const double sum = sharder.ParallelSumOverShards(
219 [&diag](
const Shard& shard) {
return shard(diag).diagonal().sum(); });
223TEST(ParallelSumOverShards, SmallExampleUsingDiagonalMatrixMultiplication) {
224 DiagonalMatrix<double, Eigen::Dynamic> diag{{1, 2, 3}};
225 VectorXd vec{{1, 1, 1}};
227 Sharder sharder(diag.cols(), 2,
nullptr);
228 sharder.ParallelForEachShard(
229 [&](
const Shard& shard) { shard(answer) = shard(diag) * shard(vec); });
230 EXPECT_THAT(answer, ElementsAre(1.0, 2.0, 3.0));
233TEST(ParallelTrueForAllShards, SmallTrueExample) {
236 Sharder sharder(vec.size(), 2,
nullptr);
237 const bool result = sharder.ParallelTrueForAllShards(
238 [&vec](
const Shard& shard) {
return (shard(vec).array() > 0.0).all(); });
242TEST(ParallelTrueForAllShards, SmallFalseExample) {
245 Sharder sharder(vec.size(), 2,
nullptr);
246 const bool result = sharder.ParallelTrueForAllShards(
247 [&vec](
const Shard& shard) {
return (shard(vec).array() < 2.5).all(); });
248 EXPECT_FALSE(result);
251TEST(MatrixVectorProductTest, SmallExample) {
252 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat =
254 Sharder sharder(mat, 3,
nullptr);
258 EXPECT_THAT(ans, ElementsAre(6.0, -0.5, 6.0, 19));
261TEST(AddScaledVectorTest, SmallExample) {
262 Sharder sharder(3, 2,
nullptr);
263 VectorXd vec1(3), vec2(3);
267 EXPECT_THAT(vec1, ElementsAre(6, 19, 26));
270TEST(AssignVectorTest, SmallExample) {
271 Sharder sharder(3, 2,
nullptr);
272 VectorXd vec1, vec2(3);
275 EXPECT_THAT(vec1, ElementsAre(1, 7, 3));
278TEST(CloneVectorTest, SmallExample) {
279 Sharder sharder(3, 2,
nullptr);
282 EXPECT_THAT(
CloneVector(vec, sharder), ElementsAre(1, 7, 3));
285TEST(CoefficientWiseProductInPlaceTest, SmallExample) {
286 Sharder sharder(3, 2,
nullptr);
287 VectorXd vec1(3), vec2(3);
292 EXPECT_THAT(vec1, ElementsAre(4, 10, 60));
295TEST(CoefficientWiseQuotientInPlaceTest, SmallExample) {
296 Sharder sharder(3, 2,
nullptr);
297 VectorXd vec1(3), vec2(3);
302 EXPECT_THAT(vec1, ElementsAre(4, 3, 4));
305TEST(DotTest, SmallExample) {
306 Sharder sharder(3, 2,
nullptr);
307 VectorXd vec1(3), vec2(3);
310 double ans =
Dot(vec1, vec2, sharder);
311 EXPECT_THAT(ans, DoubleNear(4 + 10 + 18, 1.0e-13));
314TEST(LInfNormTest, SmallExample) {
315 Sharder sharder(3, 2,
nullptr);
318 double ans =
LInfNorm(vec, sharder);
322TEST(LInfNormTest, EmptyExample) {
323 Sharder sharder(0, 2,
nullptr);
325 double ans =
LInfNorm(vec, sharder);
329TEST(L1NormTest, SmallExample) {
330 Sharder sharder(3, 2,
nullptr);
333 double ans =
L1Norm(vec, sharder);
337TEST(L1NormTest, EmptyExample) {
338 Sharder sharder(0, 2,
nullptr);
340 double ans =
L1Norm(vec, sharder);
344TEST(SquaredNormTest, SmallExample) {
345 Sharder sharder(3, 2,
nullptr);
349 EXPECT_THAT(ans, DoubleNear(1 + 4 + 9, 1.0e-13));
352TEST(NormTest, SmallExample) {
353 Sharder sharder(3, 2,
nullptr);
356 double ans =
Norm(vec, sharder);
357 EXPECT_THAT(ans, DoubleNear(std::sqrt(1 + 4 + 9), 1.0e-13));
360TEST(SquaredDistanceTest, SmallExample) {
361 Sharder sharder(3, 2,
nullptr);
367 EXPECT_THAT(ans, DoubleNear(5, 1.0e-13));
370TEST(DistanceTest, SmallExample) {
371 Sharder sharder(3, 2,
nullptr);
376 double ans =
Distance(vec1, vec2, sharder);
377 EXPECT_THAT(ans, DoubleNear(std::sqrt(5), 1.0e-13));
380TEST(ScaledLInfNormTest, SmallExample) {
381 Sharder sharder(3, 2,
nullptr);
390TEST(ScaledSquaredNormTest, SmallExample) {
391 Sharder sharder(3, 2,
nullptr);
400TEST(ScaledNormTest, SmallExample) {
401 Sharder sharder(3, 2,
nullptr);
407 EXPECT_EQ(ans, std::sqrt(169));
411 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat =
413 Sharder sharder(mat, 3,
nullptr);
414 VectorXd row_scaling_vec(3);
415 VectorXd col_scaling_vec(4);
416 row_scaling_vec << 1, -2, 1;
417 col_scaling_vec << 1, 2, -1, -1;
420 EXPECT_THAT(answer, ElementsAre(7, 1, 6, 5));
424 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat =
426 Sharder sharder(mat, 3,
nullptr);
427 VectorXd row_scaling_vec(3);
428 VectorXd col_scaling_vec(4);
429 row_scaling_vec << 1, -2, 1;
430 col_scaling_vec << 1, 2, -1, -1;
433 EXPECT_THAT(answer, ElementsAre(std::sqrt(54), 1.0, 6.0, std::sqrt(41)));
437 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat(4, 4);
438 std::vector<Eigen::Triplet<double, int64_t>> matrix_triplets = {
439 {0, 0, 1.0}, {1, 1, 2.5}, {3, 3, -3}};
440 mat.setFromTriplets(matrix_triplets.begin(), matrix_triplets.end());
441 Sharder sharder(mat, 3,
nullptr);
446 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat(3, 5);
447 std::vector<Eigen::Triplet<double, int64_t>> matrix_triplets = {
448 {0, 0, 2}, {1, 1, -1}, {2, 2, 3}};
449 mat.setFromTriplets(matrix_triplets.begin(), matrix_triplets.end());
450 Sharder sharder(mat, 3,
nullptr);
455 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat(3, 3);
456 std::vector<Eigen::Triplet<double, int64_t>> matrix_triplets = {
457 {0, 0, 2}, {0, 1, -1}, {1, 0, -1}, {2, 2, 1}};
458 mat.setFromTriplets(matrix_triplets.begin(), matrix_triplets.end());
459 Sharder sharder(mat, 3,
nullptr);
463class VariousSizesTest :
public testing::TestWithParam<int64_t> {};
465TEST_P(VariousSizesTest, LargeMatVec) {
466 const int64_t size = GetParam();
467 Eigen::SparseMatrix<double, Eigen::ColMajor, int64_t> mat =
468 LargeSparseMatrix(size);
469 const int num_threads = 5;
470 const int shards_per_thread = 3;
471 ThreadPool pool(
"MatrixVectorProductTest", num_threads);
473 Sharder sharder(mat, shards_per_thread * num_threads, &pool);
474 VectorXd rhs = VectorXd::Random(size);
475 VectorXd direct = mat.transpose() * rhs;
477 EXPECT_LE((direct - threaded).norm(), 1.0e-8);
480TEST_P(VariousSizesTest, LargeVectors) {
481 const int64_t size = GetParam();
482 const int num_threads = 5;
483 ThreadPool pool(
"SquaredNormTest", num_threads);
485 Sharder sharder(size, num_threads, &pool);
486 VectorXd vec = VectorXd::Random(size);
487 const double direct = vec.squaredNorm();
489 EXPECT_THAT(threaded, DoubleNear(direct, size * 1.0e-14));
492INSTANTIATE_TEST_SUITE_P(VariousSizesTestInstantiation, VariousSizesTest,
493 testing::Values(10, 1000, 100 * 1000));
static IntegralType CeilOfRatio(IntegralType numerator, IntegralType denominator)
double SquaredNorm(const VectorXd &vector, const Sharder &sharder)
double ScaledNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
double Dot(const VectorXd &v1, const VectorXd &v2, const Sharder &sharder)
double SquaredDistance(const VectorXd &vector1, const VectorXd &vector2, const Sharder &sharder)
double LInfNorm(const VectorXd &vector, const Sharder &sharder)
double Distance(const VectorXd &vector1, const VectorXd &vector2, const Sharder &sharder)
VectorXd TransposedMatrixVectorProduct(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &vector, const Sharder &sharder)
double ScaledLInfNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
double ScaledSquaredNorm(const VectorXd &vector, const VectorXd &scale, const Sharder &sharder)
VectorXd ScaledColLInfNorm(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &row_scaling_vec, const VectorXd &col_scaling_vec, const Sharder &sharder)
void AddScaledVector(const double scale, const VectorXd &increment, const Sharder &sharder, VectorXd &dest)
void CoefficientWiseProductInPlace(const VectorXd &scale, const Sharder &sharder, VectorXd &dest)
void CoefficientWiseQuotientInPlace(const VectorXd &scale, const Sharder &sharder, VectorXd &dest)
VectorXd ScaledColL2Norm(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const VectorXd &row_scaling_vec, const VectorXd &col_scaling_vec, const Sharder &sharder)
double L1Norm(const VectorXd &vector, const Sharder &sharder)
VectorXd CloneVector(const VectorXd &vec, const Sharder &sharder)
bool IsDiagonal(const Eigen::SparseMatrix< double, Eigen::ColMajor, int64_t > &matrix, const Sharder &sharder)
double Norm(const VectorXd &vector, const Sharder &sharder)
void AssignVector(const VectorXd &vec, const Sharder &sharder, VectorXd &dest)
TEST(LinearAssignmentTest, NullMatrix)