File: PerfTest_ViewFill.hpp

package info (click to toggle)
kokkos 4.7.01-2
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 16,636 kB
  • sloc: cpp: 223,676; sh: 2,446; makefile: 2,437; python: 91; fortran: 4; ansic: 2
file content (151 lines) | stat: -rw-r--r-- 4,180 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
//@HEADER
// ************************************************************************
//
//                        Kokkos v. 4.0
//       Copyright (2022) National Technology & Engineering
//               Solutions of Sandia, LLC (NTESS).
//
// Under the terms of Contract DE-NA0003525 with NTESS,
// the U.S. Government retains certain rights in this software.
//
// Part of Kokkos, under the Apache License v2.0 with LLVM Exceptions.
// See https://kokkos.org/LICENSE for license information.
// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
//
//@HEADER

#include "Benchmark_Context.hpp"

#include <cmath>

namespace Test {

static constexpr int N = 10;

template <class ViewType>
void fill_view(ViewType& a, typename ViewType::const_value_type& val,
               benchmark::State& state) {
  for (auto _ : state) {
    Kokkos::fence();

    Kokkos::Timer timer;
    Kokkos::deep_copy(a, val);
    KokkosBenchmark::report_results(state, a, 1, timer.seconds());
  }
}

template <
    class Layout,
    class MemorySpace = typename Kokkos::DefaultExecutionSpace::memory_space>
static void ViewFill_Rank1(benchmark::State& state) {
  const int N1 = state.range(0);
  const int N2 = N1 * N1;
  const int N4 = N2 * N2;
  const int N8 = N4 * N4;

  Kokkos::View<double*, Layout, MemorySpace> a("A1", N8);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Rank2(benchmark::State& state) {
  const int N1 = state.range(0);
  const int N2 = N1 * N1;
  const int N4 = N2 * N2;

  Kokkos::View<double**, Layout> a("A2", N4, N4);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Rank3(benchmark::State& state) {
  const int N1 = state.range(0);
  const int N2 = N1 * N1;
  const int N3 = N2 * N1;

  Kokkos::View<double***, Layout> a("A3", N3, N3, N2);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Rank4(benchmark::State& state) {
  const int N1 = state.range(0);
  const int N2 = N1 * N1;

  Kokkos::View<double****, Layout> a("A4", N2, N2, N2, N2);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Rank5(benchmark::State& state) {
  const int N1 = state.range(0);
  const int N2 = N1 * N1;

  Kokkos::View<double*****, Layout> a("A5", N2, N2, N1, N1, N2);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Rank6(benchmark::State& state) {
  const int N1 = state.range(0);
  const int N2 = N1 * N1;

  Kokkos::View<double******, Layout> a("A6", N2, N1, N1, N1, N1, N2);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Rank7(benchmark::State& state) {
  const int N1 = state.range(0);
  const int N2 = N1 * N1;

  Kokkos::View<double*******, Layout> a("A7", N2, N1, N1, N1, N1, N1, N1);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Rank8(benchmark::State& state) {
  const int N1 = state.range(0);

  Kokkos::View<double********, Layout> a("A8", N1, N1, N1, N1, N1, N1, N1, N1);
  fill_view(a, 1.1, state);
}

template <class Layout>
static void ViewFill_Raw(benchmark::State& state) {
  const int N8 = std::pow(state.range(0), 8);

  Kokkos::View<double*, Layout> a("A1", N8);
  double* a_ptr = a.data();

  for (auto _ : state) {
    Kokkos::Timer timer;
    Kokkos::parallel_for(
        N8, KOKKOS_LAMBDA(const int& i) { a_ptr[i] = 1.1; });
    Kokkos::fence();

    KokkosBenchmark::report_results(state, a, 1, timer.seconds());
  }
}

[[maybe_unused]] static void ViewFill_Rank1Strided(benchmark::State& state) {
  const size_t N8 = std::pow(state.range(0), 8);

  // This benchmark allocates more data in order to measure a view fill
  // of the same size as the contiguous benchmarks, so in cases where they
  // can be run, this one may fail to allocate data (e.g., on a small CI runner)
  try {
    // allocate 2x the size since layout only has 1/2 the elements
    Kokkos::View<double*> a("A1", N8 * 2);

    Kokkos::LayoutStride layout(N8 / 2, 2);
    Kokkos::View<double*, Kokkos::LayoutStride> a_stride(a.data(), layout);

    fill_view(a_stride, 1.1, state);

  } catch (const std::runtime_error& e) {
    state.SkipWithError(e.what());
  }
}

}  // namespace Test