Skip to content

Commit

Permalink
Benchmark update
Browse files Browse the repository at this point in the history
  • Loading branch information
mlund committed Feb 6, 2020
1 parent 6cd7138 commit 667bfab
Show file tree
Hide file tree
Showing 2 changed files with 12 additions and 16 deletions.
18 changes: 7 additions & 11 deletions src/scatter.h
Original file line number Diff line number Diff line change
Expand Up @@ -252,7 +252,7 @@ template <typename T> class SamplingPolicy {
*
* For more information, see @see http://doi.org/d8zgw5 and @see http://doi.org/10.1063/1.449987.
*/
template <typename T = double, Algorithm method = EIGEN, typename TSamplingPolicy = SamplingPolicy<T>>
template <typename T = double, Algorithm method = SIMD, typename TSamplingPolicy = SamplingPolicy<T>>
class StructureFactorPBC : private TSamplingPolicy {
//! sample directions (h,k,l)
const std::vector<Point> directions = {
Expand All @@ -276,8 +276,7 @@ class StructureFactorPBC : private TSamplingPolicy {
const Point q = (2 * pc::pi * p / boxlength) * directions[i]; // scattering vector
T sum_sin = 0.0;
T sum_cos = 0.0;
switch (method) {
case SIMD: {
if constexpr (method == SIMD) {
// When sine and cosine is computed in separate loops, advanced sine and cosine implementation
// utilizing SIMD instructions may be used to get at least 4 times performance boost.
// As of January 2020, only GCC exploits this using libmvec library if --ffast-math is enabled.
Expand All @@ -292,24 +291,21 @@ class StructureFactorPBC : private TSamplingPolicy {
for (auto &qr : qr_std) {
sum_cos += std::cos(qr);
}
} break;
case EIGEN: {
} else if constexpr (method == EIGEN) {
// Map is a Nx3 matrix facade into original std::vector
auto qr = Eigen::Map<Eigen::MatrixXd, 0, Eigen::Stride<1, 3>>((T *)positions.data(),
auto qr = Eigen::Map<Eigen::MatrixXd, 0, Eigen::Stride<1, 3>>((double *)positions.data(),
positions.size(), 3) *
q;
sum_sin = qr.array().sin().sum();
sum_cos = qr.array().cos().sum();
} break;
case GENERIC:
sum_sin = qr.array().cast<T>().sin().sum();
sum_cos = qr.array().cast<T>().cos().sum();
} else if constexpr (method == GENERIC) {
// TODO: Optimize also for other compilers than GCC by using a vector math library, e.g.,
// TODO: https://github.com/vectorclass/version2
for (auto &r : positions) { // loop over positions
T qr = q.dot(r); // scalar product q*r
sum_sin += sin(qr);
sum_cos += cos(qr);
}
break;
};
// collect average, `norm()` gives the scattering vector length
const T sf = (sum_sin * sum_sin + sum_cos * sum_cos) / (T)(positions.size());
Expand Down
10 changes: 5 additions & 5 deletions src/scatter_test.h
Original file line number Diff line number Diff line change
Expand Up @@ -32,12 +32,12 @@ TEST_CASE_TEMPLATE("[Faunus] StructureFactorPBC", T, StructureFactorPBC<double,
TEST_CASE("Benchmark") {
std::vector<Point> pos(1000);
for (auto &p : pos)
p = Eigen::Vector3d();
p = Eigen::Vector3d::Random() * box;
ankerl::nanobench::Config bench;
bench.minEpochIterations(3103);
bench.run("SIMD", [&] { StructureFactorPBC<double, SIMD>(2).sample(pos, box); }).doNotOptimizeAway();
bench.run("EIGEN", [&] { StructureFactorPBC<double, EIGEN>(2).sample(pos, box); }).doNotOptimizeAway();
bench.run("GENERIC", [&] { StructureFactorPBC<double, GENERIC>(2).sample(pos, box); }).doNotOptimizeAway();
bench.minEpochIterations(100);
bench.run("SIMD", [&] { StructureFactorPBC<double, SIMD>(10).sample(pos, box); }).doNotOptimizeAway();
bench.run("EIGEN", [&] { StructureFactorPBC<double, EIGEN>(10).sample(pos, box); }).doNotOptimizeAway();
bench.run("GENERIC", [&] { StructureFactorPBC<double, GENERIC>(10).sample(pos, box); }).doNotOptimizeAway();
}
#endif

Expand Down

0 comments on commit 667bfab

Please sign in to comment.