Page MenuHomec4science

bm_maths.cpp
No OneTemporary

File Metadata

Created
Sat, Nov 16, 02:07

bm_maths.cpp

// this file is a micro-benchmark of the maths functions
//
// It uses the google/benchmark framework
#include <benchmark/benchmark.h>
#include "specmicp_common/types.hpp"
#include "specmicp_common/eigen/incl_eigen_core.hpp"
#include <Eigen/QR>
#include <Eigen/LU>
#include "specmicp_common/eigen/incl_eigen_sparse_core.hpp"
#include <Eigen/SparseLU>
#include <Eigen/SparseQR>
#include <iostream>
#include <fstream>
static void pow_to_2(benchmark::State& state)
{
const int size {10};
Eigen::VectorXd x = Eigen::VectorXd::Random(size);
while (state.KeepRunning())
{
for (int i=0; i<size; ++i)
{
auto y = std::pow(x(i), 2);
benchmark::DoNotOptimize(y);
}
}
}
BENCHMARK(pow_to_2);
static void pow_to_3(benchmark::State& state)
{
const int size {10};
Eigen::VectorXd x = Eigen::VectorXd::Random(size);
while (state.KeepRunning())
{
for (int i=0; i<size; ++i)
{
auto y = std::pow(x(i), 3);
benchmark::DoNotOptimize(y);
}
}
}
BENCHMARK(pow_to_3);
static void pow_to_2_5(benchmark::State& state)
{
const int size {10};
Eigen::VectorXd x = Eigen::VectorXd::Random(size);
while (state.KeepRunning())
{
for (int i=0; i<size; ++i)
{
auto y = std::pow(x(i), 2.5);
benchmark::DoNotOptimize(y);
}
}
}
BENCHMARK(pow_to_2_5);
static void pow_10(benchmark::State& state)
{
const int size {10};
Eigen::VectorXd x = Eigen::VectorXd::Random(size);
while (state.KeepRunning())
{
for (int i=0; i<size; ++i)
{
auto y = std::pow(10, x(i));
benchmark::DoNotOptimize(y);
}
}
}
BENCHMARK(pow_10);
static void sqrt_2(benchmark::State& state)
{
const int size {10};
Eigen::VectorXd x = Eigen::VectorXd::Random(size);
while (state.KeepRunning())
{
for (int i=0; i<size; ++i)
{
auto y = std::sqrt(x(i));
benchmark::DoNotOptimize(y);
}
}
}
BENCHMARK(sqrt_2);
static void init_vector_10(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::VectorXd x = Eigen::VectorXd::Zero(10);
}
}
BENCHMARK(init_vector_10);
static void init_vector_50(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::VectorXd x = Eigen::VectorXd::Zero(50);
}
}
BENCHMARK(init_vector_50);
static void init_1_vector_50(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::VectorXd x = Eigen::VectorXd::Ones(50);
}
}
BENCHMARK(init_1_vector_50);
static void random_matrix_20_20(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::MatrixXd x = Eigen::MatrixXd::Random(20, 20);
benchmark::DoNotOptimize(x);
}
}
BENCHMARK(random_matrix_20_20);
static void random_vector_20(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::VectorXd x = Eigen::VectorXd::Random(20);
benchmark::DoNotOptimize(x);
}
}
BENCHMARK(random_vector_20);
static void solve_system_LU_20(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::MatrixXd A = Eigen::MatrixXd::Random(20, 20);
Eigen::VectorXd b = Eigen::VectorXd::Random(20);
Eigen::VectorXd x = A.partialPivLu().solve(b);
}
}
BENCHMARK(solve_system_LU_20);
static void solve_system_QR_20(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::MatrixXd A = Eigen::MatrixXd::Random(20, 20);
Eigen::VectorXd b = Eigen::VectorXd::Random(20);
Eigen::VectorXd x = A.colPivHouseholderQr().solve(b);
}
}
BENCHMARK(solve_system_QR_20);
static void random_matrix_50_50(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::MatrixXd x = Eigen::MatrixXd::Random(50, 50);
benchmark::DoNotOptimize(x);
}
}
BENCHMARK(random_matrix_50_50);
static void random_vector_50(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::VectorXd x = Eigen::VectorXd::Random(50);
benchmark::DoNotOptimize(x);
}
}
BENCHMARK(random_vector_50);
static void solve_system_LU_50(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::MatrixXd A = Eigen::MatrixXd::Random(50, 50);
Eigen::VectorXd b = Eigen::VectorXd::Random(50);
Eigen::VectorXd x = A.partialPivLu().solve(b);
}
}
BENCHMARK(solve_system_LU_50);
static void solve_system_QR_50(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::MatrixXd A = Eigen::MatrixXd::Random(50, 50);
Eigen::VectorXd b = Eigen::VectorXd::Random(50);
Eigen::VectorXd x = A.householderQr().solve(b);
}
}
BENCHMARK(solve_system_QR_50);
static void solve_system_colQR_50(benchmark::State& state)
{
while(state.KeepRunning())
{
Eigen::MatrixXd A = Eigen::MatrixXd::Random(50, 50);
Eigen::VectorXd b = Eigen::VectorXd::Random(50);
Eigen::VectorXd x = A.colPivHouseholderQr().solve(b);
}
}
BENCHMARK(solve_system_colQR_50);
using T = Eigen::Triplet<double>;
using SparseT = Eigen::SparseMatrix<double>;
static std::vector<T> get_triplets_list()
{
const int size = 50;
std::vector<T> list_triplets;
list_triplets.reserve(3*size);
list_triplets.push_back(T(0, 0, 20));
list_triplets.push_back(T(0, 1, -10));
for (int i=1; i<size-1; ++i)
{
list_triplets.push_back(T(i, i-1, -10));
list_triplets.push_back(T(i, i, 10));
list_triplets.push_back(T(i, i, 10));
list_triplets.push_back(T(i, i+1, -10));
}
list_triplets.push_back(T(size-1, size-1, 10));
list_triplets.push_back(T(size-1, size-2, -10));
return list_triplets;
}
static void set_from_triplets(benchmark::State& state)
{
auto list_triplets = get_triplets_list();
while(state.KeepRunning())
{
SparseT mat(50, 50);
mat.setFromTriplets(list_triplets.begin(), list_triplets.end());
}
}
BENCHMARK(set_from_triplets);
static void solve_sparse_lu(benchmark::State &state)
{
SparseT mat(50, 50);
auto list_triplets = get_triplets_list();
mat.setFromTriplets(list_triplets.begin(), list_triplets.end());
while(state.KeepRunning())
{
Eigen::VectorXd x = Eigen::VectorXd::Random(50);
Eigen::SparseLU<SparseT, Eigen::COLAMDOrdering<int>> solver;
solver.compute(mat);
Eigen::VectorXd sol = solver.solve(x);
}
}
BENCHMARK(solve_sparse_lu);
static void solve_sparse_qr(benchmark::State &state)
{
SparseT mat(50, 50);
auto list_triplets = get_triplets_list();
mat.setFromTriplets(list_triplets.begin(), list_triplets.end());
while(state.KeepRunning())
{
Eigen::VectorXd x = Eigen::VectorXd::Random(50);
Eigen::SparseQR<SparseT, Eigen::COLAMDOrdering<int>> solver;
solver.compute(mat);
Eigen::VectorXd sol = solver.solve(x);
}
}
BENCHMARK(solve_sparse_qr);
BENCHMARK_MAIN();

Event Timeline