diff --git a/examples/basic/basic.cu b/examples/basic/basic.cu index f7cc13a..4925b99 100644 --- a/examples/basic/basic.cu +++ b/examples/basic/basic.cu @@ -1,44 +1,91 @@ -#include -#include +#include "../common.h" -#include +#include #include +#include + +#include + +using cu::tangent; + +constexpr auto f(auto x, auto y) +{ + auto print = [](auto x) { printf("{%g, %g}\n", x.v, x.d); }; + + auto a = x + y; + auto b = x - y; + auto c = x * y; + auto d = x / y; + auto e = max(x, y); + auto f = min(x, y); + auto g = mid(x, y, y); + auto h = sin(x); + auto i = cos(x); + auto j = exp(x); + auto k = log(x); + auto l = pown(x, 2); -template -using mc = cu::mccormick; + print(a); + print(b); + print(c); + print(d); + print(e); + print(f); + print(g); + print(h); + print(i); + print(j); + print(k); + print(l); + return a; +} + +__global__ void kernel(tangent *xs, tangent *ys, + tangent *res, int n) +{ + int i = blockIdx.x * blockDim.x + threadIdx.x; + if (i < n) { + res[i] = f(xs[i], ys[i]); + } +} int main() { - // constexpr int n = 256; - // using T = mc; - // T xs[n], ys[n], res[n]; - // - // // generate dummy data - // for (int i = 0; i < n; i++) { - // double v = i; - // xs[i] = { .cv = -v, .cc = v, .box = { .lb = -v, .ub = v } }; - // ys[i] = { .cv = -v, .cc = v, .box = { .lb = -v, .ub = v } }; + constexpr int n = 16; + using T = tangent; + T xs[n], ys[n], res[n]; + + // generate dummy data + for (int i = 0; i < n; i++) { + double v = i; + xs[i] = { v, 1.0 }; + ys[i] = { v, 0.0 }; + } + + // for (auto el : xs) { + // std::cout << el << std::endl; + // } + + T *d_xs, *d_ys, *d_res; + CUDA_CHECK(cudaMalloc(&d_xs, n * sizeof(*xs))); + CUDA_CHECK(cudaMalloc(&d_ys, n * sizeof(*ys))); + CUDA_CHECK(cudaMalloc(&d_res, n * sizeof(*res))); + + CUDA_CHECK(cudaMemcpy(d_xs, xs, n * sizeof(*xs), cudaMemcpyHostToDevice)); + CUDA_CHECK(cudaMemcpy(d_ys, ys, n * sizeof(*ys), cudaMemcpyHostToDevice)); + + kernel<<>>(d_xs, d_ys, d_res, n); + + CUDA_CHECK(cudaMemcpy(res, d_res, n * sizeof(*res), cudaMemcpyDeviceToHost)); + + // for (auto el : res) { + // std::cout << el << std::endl; // } - // - // mc *d_xs, *d_ys, *d_res; - // CUDA_CHECK(cudaMalloc(&d_xs, n * sizeof(*xs))); - // CUDA_CHECK(cudaMalloc(&d_ys, n * sizeof(*ys))); - // CUDA_CHECK(cudaMalloc(&d_res, n * sizeof(*res))); - // - // CUDA_CHECK(cudaMemcpy(d_xs, xs, n * sizeof(*xs), cudaMemcpyHostToDevice)); - // CUDA_CHECK(cudaMemcpy(d_ys, ys, n * sizeof(*ys), cudaMemcpyHostToDevice)); - // - // kernel<<>>(d_xs, d_ys, d_res, n); - // - // CUDA_CHECK(cudaMemcpy(res, d_res, n * sizeof(*res), cudaMemcpyDeviceToHost)); - // - // auto r = res[0]; - // printf("beale(0, 0) = " MCCORMICK_FORMAT "\n", r.box.lb, r.cv, r.cc, r.box.ub); - // - // CUDA_CHECK(cudaFree(d_xs)); - // CUDA_CHECK(cudaFree(d_ys)); - // CUDA_CHECK(cudaFree(d_res)); + + CUDA_CHECK(cudaFree(d_xs)); + CUDA_CHECK(cudaFree(d_ys)); + CUDA_CHECK(cudaFree(d_res)); return 0; } diff --git a/examples/common.h b/examples/common.h new file mode 100644 index 0000000..2f23d04 --- /dev/null +++ b/examples/common.h @@ -0,0 +1,21 @@ +#ifndef CUTANGENT_COMMON_H +#define CUTANGENT_COMMON_H + +#include +#include +#include + +#include + +#define CUDA_CHECK(x) \ + do { \ + cudaError_t err = x; \ + if (err != cudaSuccess) { \ + fprintf(stderr, "CUDA error in %s at %s:%d: %s (%s=%d)\n", __FUNCTION__, \ + __FILE__, __LINE__, cudaGetErrorString(err), \ + cudaGetErrorName(err), err); \ + abort(); \ + } \ + } while (0) + +#endif // CUTANGENT_COMMON_H diff --git a/include/cutangent/arithmetic/basic.cuh b/include/cutangent/arithmetic/basic.cuh index 60c3ad9..7c0315e 100644 --- a/include/cutangent/arithmetic/basic.cuh +++ b/include/cutangent/arithmetic/basic.cuh @@ -4,28 +4,45 @@ #include #include -// #include -// #include +#include namespace cu { -#define cuda_fn inline constexpr __device__ +#define fn inline constexpr __device__ template -cuda_fn tangent operator+(tangent a, tangent b) +fn tangent operator-(tangent x) +{ + return { -x.v, -x.d }; +} + +template +fn tangent operator+(tangent a, tangent b) { return { .v = a.v + b.v, .d = a.d + b.d }; } template -cuda_fn tangent operator*(tangent a, tangent b) +fn tangent operator-(tangent a, tangent b) +{ + return { .v = a.v - b.v, .d = a.d - b.d }; +} + +template +fn tangent operator*(tangent a, tangent b) { return { .v = a.v * b.v, .d = a.v * b.d + a.d * b.v }; } template -cuda_fn tangent max(tangent a, tangent b) +fn tangent operator/(tangent a, tangent b) +{ + return { .v = a.v / b.v, .d = (a.d * b.v - a.v * b.d) / (b.v * b.v) }; +} + +template +fn tangent max(tangent a, tangent b) { using std::max; @@ -34,7 +51,7 @@ cuda_fn tangent max(tangent a, tangent b) } template -cuda_fn tangent min(tangent a, tangent b) +fn tangent min(tangent a, tangent b) { using std::min; @@ -43,11 +60,69 @@ cuda_fn tangent min(tangent a, tangent b) } template -cuda_fn tangent mid(tangent v, tangent lb, tangent ub) +fn tangent clamp(tangent v, tangent lb, tangent ub) { - return { .v = mid(v.v, lb.v, ub.v), .d = lb.d * (v.v < lb.v) + ub.d * (v.v > ub.v) }; + using std::clamp; + + return { .v = clamp(v.v, lb.v, ub.v), .d = lb.d * (v.v < lb.v) + ub.d * (v.v > ub.v) }; +} + +template +fn tangent mid(tangent v, tangent lb, tangent ub) +{ + return clamp(v, lb, ub); +} + +template +fn tangent sin(tangent x) +{ + using std::cos; + using std::sin; + + return { .v = sin(x.v), .d = cos(x.v) * x.d }; } +template +fn tangent cos(tangent x) +{ + using std::cos; + using std::sin; + + return { .v = cos(x.v), .d = -sin(x.v) * x.d }; +} + +template +fn tangent exp(tangent x) +{ + using std::exp; + + return { .v = exp(x.v), .d = exp(x.v) * x.d }; +} + +template +fn tangent log(tangent x) +{ + using std::log; + + return { .v = log(x.v), .d = x.d / x.v }; +} + +template +fn tangent sqr(tangent x) +{ + return { .v = sqr(x.v), .d = 2.0 * x.v * x.d }; +} + +template +fn tangent pown(tangent x, auto n) +{ + using std::pow; + + return { .v = pow(x.v, n), .d = n * pow(x.v, n - 1) * x.d }; +} + +#undef fn + } // namespace cu #endif // CUTANGENT_ARITHMETIC_BASIC_CUH diff --git a/include/cutangent/format.h b/include/cutangent/format.h new file mode 100644 index 0000000..b1100e9 --- /dev/null +++ b/include/cutangent/format.h @@ -0,0 +1,19 @@ +#ifndef CUTANGENT_FORMAT_H +#define CUTANGENT_FORMAT_H + +#include + +#include + +namespace cu +{ + +template +std::ostream &operator<<(std::ostream &os, tangent x) +{ + return os << "{v: " << x.v << ", d: " << x.d << "}"; +} + +} // namespace cu + +#endif // CUTANGENT_FORMAT_H diff --git a/include/cutangent/tangent.h b/include/cutangent/tangent.h index 4ade649..4b02c56 100644 --- a/include/cutangent/tangent.h +++ b/include/cutangent/tangent.h @@ -9,12 +9,13 @@ namespace cu template struct tangent { - T v; - T d; + T v; // value + T d; // derivative constexpr auto operator<=>(const tangent &) const = default; }; + } // namespace cu #endif // CUTANGENT_TANGENT_H