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 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210
|
/*
* GridTools
*
* Copyright (c) 2014-2023, ETH Zurich
* All rights reserved.
*
* Please, refer to the LICENSE file in the root directory.
* SPDX-License-Identifier: BSD-3-Clause
*/
#include "gtest/gtest.h"
#include <cstdlib>
#include <gridtools/common/atomic_functions.hpp>
#include <gridtools/common/cuda_util.hpp>
#include <gridtools/common/defs.hpp>
template <typename T>
struct verifier {
static void TestEQ(T val, T exp) {
T err = std::fabs(val - exp) / std::fabs(val);
ASSERT_TRUE(err < 1e-12);
}
};
template <>
struct verifier<float> {
static void TestEQ(float val, float exp) {
double err = std::fabs(val - exp) / std::fabs(val);
ASSERT_TRUE(err < 1e-6);
}
};
template <>
struct verifier<int> {
static void TestEQ(int val, int exp) { ASSERT_EQ(val, exp); }
};
template <typename T>
__global__ void atomic_add_kernel(T *pReduced, const T *field, const int size) {
const int i = threadIdx.x + blockIdx.x * blockDim.x;
const int j = threadIdx.y + blockIdx.y * blockDim.y;
const int pos = j * gridDim.x * blockDim.x + i;
gridtools::atomic_add(*pReduced, field[pos]);
}
template <typename T>
__global__ void atomic_sub_kernel(T *pReduced, const T *field, const int size) {
const int i = threadIdx.x + blockIdx.x * blockDim.x;
const int j = threadIdx.y + blockIdx.y * blockDim.y;
const int pos = j * gridDim.x * blockDim.x + i;
gridtools::atomic_sub(*pReduced, field[pos]);
}
template <typename T>
__global__ void atomic_min_kernel(T *pReduced, const T *field, const int size) {
const int i = threadIdx.x + blockIdx.x * blockDim.x;
const int j = threadIdx.y + blockIdx.y * blockDim.y;
const int pos = j * gridDim.x * blockDim.x + i;
gridtools::atomic_min(*pReduced, field[pos]);
}
template <typename T>
__global__ void atomic_max_kernel(T *pReduced, const T *field, const int size) {
const int i = threadIdx.x + blockIdx.x * blockDim.x;
const int j = threadIdx.y + blockIdx.y * blockDim.y;
const int pos = j * gridDim.x * blockDim.x + i;
gridtools::atomic_max(*pReduced, field[pos]);
}
template <typename T>
void test_atomic_add() {
dim3 threadsPerBlock(4, 4);
dim3 numberOfBlocks(4, 4);
int size = threadsPerBlock.x * threadsPerBlock.y * numberOfBlocks.x * numberOfBlocks.y;
T field[size];
T sumRef = 0;
T sum = 0;
T *sumDevice;
GT_CUDA_CHECK(cudaMalloc(&sumDevice, sizeof(T)));
GT_CUDA_CHECK(cudaMemcpy(sumDevice, &sum, sizeof(T), cudaMemcpyHostToDevice));
T *fieldDevice;
GT_CUDA_CHECK(cudaMalloc(&fieldDevice, sizeof(T) * size));
for (int cnt = 0; cnt < size; ++cnt) {
field[cnt] = static_cast<T>(std::rand() % 100 + (std::rand() % 100) * 0.005);
sumRef += field[cnt];
}
GT_CUDA_CHECK(cudaMemcpy(fieldDevice, &field[0], sizeof(T) * size, cudaMemcpyHostToDevice));
atomic_add_kernel<<<numberOfBlocks, threadsPerBlock>>>(sumDevice, fieldDevice, size);
GT_CUDA_CHECK(cudaMemcpy(&sum, sumDevice, sizeof(T), cudaMemcpyDeviceToHost));
verifier<T>::TestEQ(sumRef, sum);
}
template <typename T>
void test_atomic_sub() {
dim3 threadsPerBlock(4, 4);
dim3 numberOfBlocks(4, 4);
int size = threadsPerBlock.x * threadsPerBlock.y * numberOfBlocks.x * numberOfBlocks.y;
T field[size];
T sumRef = 0;
T sum = 0;
T *sumDevice;
GT_CUDA_CHECK(cudaMalloc(&sumDevice, sizeof(T)));
GT_CUDA_CHECK(cudaMemcpy(sumDevice, &sum, sizeof(T), cudaMemcpyHostToDevice));
T *fieldDevice;
GT_CUDA_CHECK(cudaMalloc(&fieldDevice, sizeof(T) * size));
for (int cnt = 0; cnt < size; ++cnt) {
field[cnt] = static_cast<T>(std::rand() % 100 + (std::rand() % 100) * 0.005);
sumRef -= field[cnt];
}
GT_CUDA_CHECK(cudaMemcpy(fieldDevice, &field[0], sizeof(T) * size, cudaMemcpyHostToDevice));
atomic_sub_kernel<<<numberOfBlocks, threadsPerBlock>>>(sumDevice, fieldDevice, size);
GT_CUDA_CHECK(cudaMemcpy(&sum, sumDevice, sizeof(T), cudaMemcpyDeviceToHost));
verifier<T>::TestEQ(sumRef, sum);
}
template <typename T>
void test_atomic_min() {
dim3 threadsPerBlock(4, 4);
dim3 numberOfBlocks(4, 4);
int size = threadsPerBlock.x * threadsPerBlock.y * numberOfBlocks.x * numberOfBlocks.y;
T field[size];
T minRef = 99999;
T min = 99999;
T *minDevice;
GT_CUDA_CHECK(cudaMalloc(&minDevice, sizeof(T)));
GT_CUDA_CHECK(cudaMemcpy(minDevice, &min, sizeof(T), cudaMemcpyHostToDevice));
T *fieldDevice;
GT_CUDA_CHECK(cudaMalloc(&fieldDevice, sizeof(T) * size));
for (int cnt = 0; cnt < size; ++cnt) {
field[cnt] = static_cast<T>(std::rand() % 100 + (std::rand() % 100) * 0.005);
minRef = std::min(minRef, field[cnt]);
}
GT_CUDA_CHECK(cudaMemcpy(fieldDevice, &field[0], sizeof(T) * size, cudaMemcpyHostToDevice));
atomic_min_kernel<<<numberOfBlocks, threadsPerBlock>>>(minDevice, fieldDevice, size);
GT_CUDA_CHECK(cudaMemcpy(&min, minDevice, sizeof(T), cudaMemcpyDeviceToHost));
verifier<T>::TestEQ(minRef, min);
}
template <typename T>
void test_atomic_max() {
dim3 threadsPerBlock(4, 4);
dim3 numberOfBlocks(4, 4);
int size = threadsPerBlock.x * threadsPerBlock.y * numberOfBlocks.x * numberOfBlocks.y;
T field[size];
T maxRef = -1;
T max = -1;
T *maxDevice;
GT_CUDA_CHECK(cudaMalloc(&maxDevice, sizeof(T)));
GT_CUDA_CHECK(cudaMemcpy(maxDevice, &max, sizeof(T), cudaMemcpyHostToDevice));
T *fieldDevice;
GT_CUDA_CHECK(cudaMalloc(&fieldDevice, sizeof(T) * size));
for (int cnt = 0; cnt < size; ++cnt) {
field[cnt] = static_cast<T>(std::rand() % 100 + (std::rand() % 100) * 0.005);
maxRef = std::max(maxRef, field[cnt]);
}
GT_CUDA_CHECK(cudaMemcpy(fieldDevice, &field[0], sizeof(T) * size, cudaMemcpyHostToDevice));
atomic_max_kernel<<<numberOfBlocks, threadsPerBlock>>>(maxDevice, fieldDevice, size);
GT_CUDA_CHECK(cudaMemcpy(&max, maxDevice, sizeof(T), cudaMemcpyDeviceToHost));
verifier<T>::TestEQ(maxRef, max);
}
TEST(AtomicFunctionsUnittest, atomic_add_int) { test_atomic_add<int>(); }
TEST(AtomicFunctionsUnittest, atomic_add_real) {
test_atomic_add<double>();
test_atomic_add<float>();
}
TEST(AtomicFunctionsUnittest, atomic_sub_int) { test_atomic_sub<int>(); }
TEST(AtomicFunctionsUnittest, atomic_sub_real) {
test_atomic_sub<double>();
test_atomic_sub<float>();
}
TEST(AtomicFunctionsUnittest, atomic_min_int) { test_atomic_min<int>(); }
TEST(AtomicFunctionsUnittest, atomic_min_real) {
test_atomic_min<double>();
test_atomic_min<float>();
}
TEST(AtomicFunctionsUnittest, atomic_max_int) { test_atomic_max<int>(); }
TEST(AtomicFunctionsUnittest, atomic_max_real) {
test_atomic_max<double>();
test_atomic_max<float>();
}
|