xref: /libCEED/include/ceed/jit-source/cuda/cuda-ref-basis-tensor.h (revision 0b2e49137516da01e333647ed4c60bbf15a6a3be)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2a0154adeSJed Brown // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a0154adeSJed Brown //
4a0154adeSJed Brown // SPDX-License-Identifier: BSD-2-Clause
5a0154adeSJed Brown //
6a0154adeSJed Brown // This file is part of CEED:  http://github.com/ceed
7a0154adeSJed Brown 
8b2165e7aSSebastian Grimberg /// @file
9b2165e7aSSebastian Grimberg /// Internal header for CUDA tensor product basis
10b2165e7aSSebastian Grimberg 
11c9c2c079SJeremy L Thompson #include <ceed.h>
12a0154adeSJed Brown 
13a0154adeSJed Brown //------------------------------------------------------------------------------
14a0154adeSJed Brown // Tensor Basis Kernels
15a0154adeSJed Brown //------------------------------------------------------------------------------
16a0154adeSJed Brown 
17a0154adeSJed Brown //------------------------------------------------------------------------------
18a0154adeSJed Brown // Interp
19a0154adeSJed Brown //------------------------------------------------------------------------------
20*0b2e4913SJeremy L Thompson extern "C" __global__ void Interp(const CeedInt num_elem, const CeedInt is_transpose, const CeedScalar *__restrict__ interp_1d,
212b730f8bSJeremy L Thompson                                   const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
22a0154adeSJed Brown   const CeedInt i = threadIdx.x;
23a0154adeSJed Brown 
24a0154adeSJed Brown   __shared__ CeedScalar s_mem[BASIS_Q_1D * BASIS_P_1D + 2 * BASIS_BUF_LEN];
25a0154adeSJed Brown   CeedScalar           *s_interp_1d = s_mem;
26a0154adeSJed Brown   CeedScalar           *s_buffer_1  = s_mem + BASIS_Q_1D * BASIS_P_1D;
27a0154adeSJed Brown   CeedScalar           *s_buffer_2  = s_buffer_1 + BASIS_BUF_LEN;
28a0154adeSJed Brown   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
29a0154adeSJed Brown     s_interp_1d[k] = interp_1d[k];
30a0154adeSJed Brown   }
31a0154adeSJed Brown 
32*0b2e4913SJeremy L Thompson   const CeedInt P             = is_transpose ? BASIS_Q_1D : BASIS_P_1D;
33*0b2e4913SJeremy L Thompson   const CeedInt Q             = is_transpose ? BASIS_P_1D : BASIS_Q_1D;
34*0b2e4913SJeremy L Thompson   const CeedInt stride_0      = is_transpose ? 1 : BASIS_P_1D;
35*0b2e4913SJeremy L Thompson   const CeedInt stride_1      = is_transpose ? BASIS_P_1D : 1;
36*0b2e4913SJeremy L Thompson   const CeedInt u_stride      = is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
37*0b2e4913SJeremy L Thompson   const CeedInt v_stride      = is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
38*0b2e4913SJeremy L Thompson   const CeedInt u_comp_stride = num_elem * (is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
39*0b2e4913SJeremy L Thompson   const CeedInt v_comp_stride = num_elem * (is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
40*0b2e4913SJeremy L Thompson   const CeedInt u_size        = is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
41a0154adeSJed Brown 
42a0154adeSJed Brown   // Apply basis element by element
43a0154adeSJed Brown   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
44a0154adeSJed Brown     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
45a0154adeSJed Brown       const CeedScalar *cur_u = u + elem * u_stride + comp * u_comp_stride;
46a0154adeSJed Brown       CeedScalar       *cur_v = v + elem * v_stride + comp * v_comp_stride;
47672b0f2aSSebastian Grimberg       CeedInt           pre   = u_size;
48672b0f2aSSebastian Grimberg       CeedInt           post  = 1;
49ca735530SJeremy L Thompson 
50a0154adeSJed Brown       for (CeedInt k = i; k < u_size; k += blockDim.x) {
51a0154adeSJed Brown         s_buffer_1[k] = cur_u[k];
52a0154adeSJed Brown       }
53a0154adeSJed Brown       for (CeedInt d = 0; d < BASIS_DIM; d++) {
54a0154adeSJed Brown         __syncthreads();
552a86cc9dSSebastian Grimberg         // Update buffers used
56a0154adeSJed Brown         pre /= P;
57a0154adeSJed Brown         const CeedScalar *in       = d % 2 ? s_buffer_2 : s_buffer_1;
58a0154adeSJed Brown         CeedScalar       *out      = d == BASIS_DIM - 1 ? cur_v : (d % 2 ? s_buffer_1 : s_buffer_2);
59672b0f2aSSebastian Grimberg         const CeedInt     writeLen = pre * post * Q;
60a0154adeSJed Brown 
61a0154adeSJed Brown         // Contract along middle index
62a0154adeSJed Brown         for (CeedInt k = i; k < writeLen; k += blockDim.x) {
63a0154adeSJed Brown           const CeedInt c  = k % post;
64a0154adeSJed Brown           const CeedInt j  = (k / post) % Q;
65a0154adeSJed Brown           const CeedInt a  = k / (post * Q);
66a0154adeSJed Brown           CeedScalar    vk = 0;
67672b0f2aSSebastian Grimberg 
682b730f8bSJeremy L Thompson           for (CeedInt b = 0; b < P; b++) vk += s_interp_1d[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c];
69a0154adeSJed Brown           out[k] = vk;
70a0154adeSJed Brown         }
71a0154adeSJed Brown         post *= Q;
72a0154adeSJed Brown       }
73a0154adeSJed Brown     }
74a0154adeSJed Brown   }
75a0154adeSJed Brown }
76a0154adeSJed Brown 
77a0154adeSJed Brown //------------------------------------------------------------------------------
78a0154adeSJed Brown // Grad
79a0154adeSJed Brown //------------------------------------------------------------------------------
80*0b2e4913SJeremy L Thompson extern "C" __global__ void Grad(const CeedInt num_elem, const CeedInt is_transpose, const CeedScalar *__restrict__ interp_1d,
812b730f8bSJeremy L Thompson                                 const CeedScalar *__restrict__ grad_1d, const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
82a0154adeSJed Brown   const CeedInt i = threadIdx.x;
83a0154adeSJed Brown 
84a0154adeSJed Brown   __shared__ CeedScalar s_mem[2 * (BASIS_Q_1D * BASIS_P_1D + BASIS_BUF_LEN)];
85a0154adeSJed Brown   CeedScalar           *s_interp_1d = s_mem;
86a0154adeSJed Brown   CeedScalar           *s_grad_1d   = s_interp_1d + BASIS_Q_1D * BASIS_P_1D;
87a0154adeSJed Brown   CeedScalar           *s_buffer_1  = s_grad_1d + BASIS_Q_1D * BASIS_P_1D;
88a0154adeSJed Brown   CeedScalar           *s_buffer_2  = s_buffer_1 + BASIS_BUF_LEN;
89a0154adeSJed Brown   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
90a0154adeSJed Brown     s_interp_1d[k] = interp_1d[k];
91a0154adeSJed Brown     s_grad_1d[k]   = grad_1d[k];
92a0154adeSJed Brown   }
93a0154adeSJed Brown 
94*0b2e4913SJeremy L Thompson   const CeedInt P             = is_transpose ? BASIS_Q_1D : BASIS_P_1D;
95*0b2e4913SJeremy L Thompson   const CeedInt Q             = is_transpose ? BASIS_P_1D : BASIS_Q_1D;
96*0b2e4913SJeremy L Thompson   const CeedInt stride_0      = is_transpose ? 1 : BASIS_P_1D;
97*0b2e4913SJeremy L Thompson   const CeedInt stride_1      = is_transpose ? BASIS_P_1D : 1;
98*0b2e4913SJeremy L Thompson   const CeedInt u_stride      = is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
99*0b2e4913SJeremy L Thompson   const CeedInt v_stride      = is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
100*0b2e4913SJeremy L Thompson   const CeedInt u_comp_stride = num_elem * (is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
101*0b2e4913SJeremy L Thompson   const CeedInt v_comp_stride = num_elem * (is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
102*0b2e4913SJeremy L Thompson   const CeedInt u_dim_stride  = is_transpose ? num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP : 0;
103*0b2e4913SJeremy L Thompson   const CeedInt v_dim_stride  = is_transpose ? 0 : num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP;
104a0154adeSJed Brown 
105a0154adeSJed Brown   // Apply basis element by element
106a0154adeSJed Brown   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
107a0154adeSJed Brown     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
108a0154adeSJed Brown       // dim*dim contractions for grad
109a0154adeSJed Brown       for (CeedInt dim_1 = 0; dim_1 < BASIS_DIM; dim_1++) {
110*0b2e4913SJeremy L Thompson         CeedInt           pre   = is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
111a0154adeSJed Brown         CeedInt           post  = 1;
1122b730f8bSJeremy L Thompson         const CeedScalar *cur_u = u + elem * u_stride + dim_1 * u_dim_stride + comp * u_comp_stride;
1132b730f8bSJeremy L Thompson         CeedScalar       *cur_v = v + elem * v_stride + dim_1 * v_dim_stride + comp * v_comp_stride;
114ca735530SJeremy L Thompson 
115a0154adeSJed Brown         for (CeedInt dim_2 = 0; dim_2 < BASIS_DIM; dim_2++) {
116a0154adeSJed Brown           __syncthreads();
117b2165e7aSSebastian Grimberg           // Update buffers used
118a0154adeSJed Brown           pre /= P;
119a0154adeSJed Brown           const CeedScalar *op       = dim_1 == dim_2 ? s_grad_1d : s_interp_1d;
1202b730f8bSJeremy L Thompson           const CeedScalar *in       = dim_2 == 0 ? cur_u : (dim_2 % 2 ? s_buffer_2 : s_buffer_1);
1212b730f8bSJeremy L Thompson           CeedScalar       *out      = dim_2 == BASIS_DIM - 1 ? cur_v : (dim_2 % 2 ? s_buffer_1 : s_buffer_2);
122672b0f2aSSebastian Grimberg           const CeedInt     writeLen = pre * post * Q;
123a0154adeSJed Brown 
124a0154adeSJed Brown           // Contract along middle index
125a0154adeSJed Brown           for (CeedInt k = i; k < writeLen; k += blockDim.x) {
126a0154adeSJed Brown             const CeedInt c   = k % post;
127a0154adeSJed Brown             const CeedInt j   = (k / post) % Q;
128a0154adeSJed Brown             const CeedInt a   = k / (post * Q);
129a0154adeSJed Brown             CeedScalar    v_k = 0;
130ca735530SJeremy L Thompson 
1312b730f8bSJeremy L Thompson             for (CeedInt b = 0; b < P; b++) v_k += op[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c];
132*0b2e4913SJeremy L Thompson             if (is_transpose && dim_2 == BASIS_DIM - 1) out[k] += v_k;
1332b730f8bSJeremy L Thompson             else out[k] = v_k;
134a0154adeSJed Brown           }
135a0154adeSJed Brown           post *= Q;
136a0154adeSJed Brown         }
137a0154adeSJed Brown       }
138a0154adeSJed Brown     }
139a0154adeSJed Brown   }
140a0154adeSJed Brown }
141a0154adeSJed Brown 
142a0154adeSJed Brown //------------------------------------------------------------------------------
143a0154adeSJed Brown // 1D quadrature weights
144a0154adeSJed Brown //------------------------------------------------------------------------------
1452b730f8bSJeremy L Thompson __device__ void Weight1d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
146a0154adeSJed Brown   const CeedInt i = threadIdx.x;
147ca735530SJeremy L Thompson 
148a0154adeSJed Brown   if (i < BASIS_Q_1D) {
149a0154adeSJed Brown     const size_t elem = blockIdx.x;
150ca735530SJeremy L Thompson 
1512b730f8bSJeremy L Thompson     if (elem < num_elem) w[elem * BASIS_Q_1D + i] = q_weight_1d[i];
152a0154adeSJed Brown   }
153a0154adeSJed Brown }
154a0154adeSJed Brown 
155a0154adeSJed Brown //------------------------------------------------------------------------------
156a0154adeSJed Brown // 2D quadrature weights
157a0154adeSJed Brown //------------------------------------------------------------------------------
1582b730f8bSJeremy L Thompson __device__ void Weight2d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
159a0154adeSJed Brown   const CeedInt i = threadIdx.x;
160a0154adeSJed Brown   const CeedInt j = threadIdx.y;
161ca735530SJeremy L Thompson 
162a0154adeSJed Brown   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
163a0154adeSJed Brown     const size_t elem = blockIdx.x;
164ca735530SJeremy L Thompson 
165a0154adeSJed Brown     if (elem < num_elem) {
166a0154adeSJed Brown       const size_t ind = (elem * BASIS_Q_1D + j) * BASIS_Q_1D + i;
167ca735530SJeremy L Thompson 
168a0154adeSJed Brown       w[ind] = q_weight_1d[i] * q_weight_1d[j];
169a0154adeSJed Brown     }
170a0154adeSJed Brown   }
171a0154adeSJed Brown }
172a0154adeSJed Brown 
173a0154adeSJed Brown //------------------------------------------------------------------------------
174a0154adeSJed Brown // 3D quadrature weights
175a0154adeSJed Brown //------------------------------------------------------------------------------
1762b730f8bSJeremy L Thompson __device__ void Weight3d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
177a0154adeSJed Brown   const CeedInt i = threadIdx.x;
178a0154adeSJed Brown   const CeedInt j = threadIdx.y;
179ca735530SJeremy L Thompson 
180a0154adeSJed Brown   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
181a0154adeSJed Brown     const size_t elem = blockIdx.x;
182ca735530SJeremy L Thompson 
183a0154adeSJed Brown     if (elem < num_elem) {
184a0154adeSJed Brown       for (CeedInt k = 0; k < BASIS_Q_1D; k++) {
185a0154adeSJed Brown         const size_t ind = ((elem * BASIS_Q_1D + k) * BASIS_Q_1D + j) * BASIS_Q_1D + i;
186ca735530SJeremy L Thompson 
187a0154adeSJed Brown         w[ind] = q_weight_1d[i] * q_weight_1d[j] * q_weight_1d[k];
188a0154adeSJed Brown       }
189a0154adeSJed Brown     }
190a0154adeSJed Brown   }
191a0154adeSJed Brown }
192a0154adeSJed Brown 
193a0154adeSJed Brown //------------------------------------------------------------------------------
194a0154adeSJed Brown // Quadrature weights
195a0154adeSJed Brown //------------------------------------------------------------------------------
1962b730f8bSJeremy L Thompson extern "C" __global__ void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ v) {
1972b730f8bSJeremy L Thompson   if (BASIS_DIM == 1) Weight1d(num_elem, q_weight_1d, v);
1982b730f8bSJeremy L Thompson   else if (BASIS_DIM == 2) Weight2d(num_elem, q_weight_1d, v);
1992b730f8bSJeremy L Thompson   else if (BASIS_DIM == 3) Weight3d(num_elem, q_weight_1d, v);
200a0154adeSJed Brown }
201