xref: /libCEED/include/ceed/jit-source/cuda/cuda-ref-basis-tensor.h (revision c9c2c07970382857cc7b4a28d359710237b91a3e)
1a0154adeSJed Brown // Copyright (c) 2017-2022, 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 
8*c9c2c079SJeremy L Thompson #include <ceed.h>
9a0154adeSJed Brown 
10a0154adeSJed Brown //------------------------------------------------------------------------------
11a0154adeSJed Brown // Tensor Basis Kernels
12a0154adeSJed Brown //------------------------------------------------------------------------------
13a0154adeSJed Brown 
14a0154adeSJed Brown //------------------------------------------------------------------------------
15a0154adeSJed Brown // Interp
16a0154adeSJed Brown //------------------------------------------------------------------------------
17a0154adeSJed Brown extern "C" __global__ void Interp(const CeedInt num_elem, const CeedInt transpose,
18a0154adeSJed Brown                                   const CeedScalar *__restrict__ interp_1d,
19a0154adeSJed Brown                                   const CeedScalar *__restrict__ u,
20a0154adeSJed Brown                                   CeedScalar *__restrict__ v) {
21a0154adeSJed Brown   const CeedInt i = threadIdx.x;
22a0154adeSJed Brown 
23a0154adeSJed Brown   __shared__ CeedScalar s_mem[BASIS_Q_1D * BASIS_P_1D + 2 * BASIS_BUF_LEN];
24a0154adeSJed Brown   CeedScalar *s_interp_1d = s_mem;
25a0154adeSJed Brown   CeedScalar *s_buffer_1 = s_mem + BASIS_Q_1D * BASIS_P_1D;
26a0154adeSJed Brown   CeedScalar *s_buffer_2 = s_buffer_1 + BASIS_BUF_LEN;
27a0154adeSJed Brown   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
28a0154adeSJed Brown     s_interp_1d[k] = interp_1d[k];
29a0154adeSJed Brown   }
30a0154adeSJed Brown 
31a0154adeSJed Brown   const CeedInt P = transpose ? BASIS_Q_1D : BASIS_P_1D;
32a0154adeSJed Brown   const CeedInt Q = transpose ? BASIS_P_1D : BASIS_Q_1D;
33a0154adeSJed Brown   const CeedInt stride_0 = transpose ? 1 : BASIS_P_1D;
34a0154adeSJed Brown   const CeedInt stride_1 = transpose ? BASIS_P_1D : 1;
35a0154adeSJed Brown   const CeedInt u_stride = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
36a0154adeSJed Brown   const CeedInt v_stride = transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
37a0154adeSJed Brown   const CeedInt u_comp_stride = num_elem * (transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
38a0154adeSJed Brown   const CeedInt v_comp_stride = num_elem * (transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
39a0154adeSJed Brown   const CeedInt u_size = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
40a0154adeSJed Brown 
41a0154adeSJed Brown   // Apply basis element by element
42a0154adeSJed Brown   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
43a0154adeSJed Brown     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
44a0154adeSJed Brown       const CeedScalar *cur_u = u + elem * u_stride + comp * u_comp_stride;
45a0154adeSJed Brown       CeedScalar *cur_v = v + elem * v_stride + comp * v_comp_stride;
46a0154adeSJed Brown       for (CeedInt k = i; k < u_size; k += blockDim.x) {
47a0154adeSJed Brown         s_buffer_1[k] = cur_u[k];
48a0154adeSJed Brown       }
49a0154adeSJed Brown       CeedInt pre = u_size;
50a0154adeSJed Brown       CeedInt post = 1;
51a0154adeSJed Brown       for (CeedInt d = 0; d < BASIS_DIM; d++) {
52a0154adeSJed Brown         __syncthreads();
53a0154adeSJed Brown         // Update bufferfers used
54a0154adeSJed Brown         pre /= P;
55a0154adeSJed Brown         const CeedScalar *in = d % 2 ? s_buffer_2 : s_buffer_1;
56a0154adeSJed Brown         CeedScalar *out = d == BASIS_DIM - 1 ? cur_v : (d % 2 ? s_buffer_1 : s_buffer_2);
57a0154adeSJed Brown 
58a0154adeSJed Brown         // Contract along middle index
59a0154adeSJed Brown         const CeedInt writeLen = pre * post * Q;
60a0154adeSJed Brown         for (CeedInt k = i; k < writeLen; k += blockDim.x) {
61a0154adeSJed Brown           const CeedInt c = k % post;
62a0154adeSJed Brown           const CeedInt j = (k / post) % Q;
63a0154adeSJed Brown           const CeedInt a = k / (post * Q);
64a0154adeSJed Brown 
65a0154adeSJed Brown           CeedScalar vk = 0;
66a0154adeSJed Brown           for (CeedInt b = 0; b < P; b++)
67a0154adeSJed Brown             vk += s_interp_1d[j*stride_0 + b*stride_1] * in[(a*P + b)*post + c];
68a0154adeSJed Brown 
69a0154adeSJed Brown           out[k] = vk;
70a0154adeSJed Brown         }
71a0154adeSJed Brown 
72a0154adeSJed Brown         post *= Q;
73a0154adeSJed Brown       }
74a0154adeSJed Brown     }
75a0154adeSJed Brown   }
76a0154adeSJed Brown }
77a0154adeSJed Brown 
78a0154adeSJed Brown //------------------------------------------------------------------------------
79a0154adeSJed Brown // Grad
80a0154adeSJed Brown //------------------------------------------------------------------------------
81a0154adeSJed Brown extern "C" __global__ void Grad(const CeedInt num_elem, const CeedInt transpose,
82a0154adeSJed Brown                                 const CeedScalar *__restrict__ interp_1d,
83a0154adeSJed Brown                                 const CeedScalar *__restrict__ grad_1d,
84a0154adeSJed Brown                                 const CeedScalar *__restrict__ u,
85a0154adeSJed Brown                                 CeedScalar *__restrict__ v) {
86a0154adeSJed Brown   const CeedInt i = threadIdx.x;
87a0154adeSJed Brown 
88a0154adeSJed Brown   __shared__ CeedScalar s_mem[2 * (BASIS_Q_1D * BASIS_P_1D + BASIS_BUF_LEN)];
89a0154adeSJed Brown   CeedScalar *s_interp_1d = s_mem;
90a0154adeSJed Brown   CeedScalar *s_grad_1d = s_interp_1d + BASIS_Q_1D * BASIS_P_1D;
91a0154adeSJed Brown   CeedScalar *s_buffer_1 = s_grad_1d + BASIS_Q_1D * BASIS_P_1D;
92a0154adeSJed Brown   CeedScalar *s_buffer_2 = s_buffer_1 + BASIS_BUF_LEN;
93a0154adeSJed Brown   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
94a0154adeSJed Brown     s_interp_1d[k] = interp_1d[k];
95a0154adeSJed Brown     s_grad_1d[k] = grad_1d[k];
96a0154adeSJed Brown   }
97a0154adeSJed Brown 
98a0154adeSJed Brown   const CeedInt P = transpose ? BASIS_Q_1D : BASIS_P_1D;
99a0154adeSJed Brown   const CeedInt Q = transpose ? BASIS_P_1D : BASIS_Q_1D;
100a0154adeSJed Brown   const CeedInt stride_0 = transpose ? 1 : BASIS_P_1D;
101a0154adeSJed Brown   const CeedInt stride_1 = transpose ? BASIS_P_1D : 1;
102a0154adeSJed Brown   const CeedInt u_stride = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
103a0154adeSJed Brown   const CeedInt v_stride = transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
104a0154adeSJed Brown   const CeedInt u_comp_stride = num_elem * (transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
105a0154adeSJed Brown   const CeedInt v_comp_stride = num_elem * (transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
106a0154adeSJed Brown   const CeedInt u_dim_stride = transpose ? num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP : 0;
107a0154adeSJed Brown   const CeedInt v_dim_stride = transpose ? 0 : num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP;
108a0154adeSJed Brown 
109a0154adeSJed Brown   // Apply basis element by element
110a0154adeSJed Brown   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
111a0154adeSJed Brown     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
112a0154adeSJed Brown 
113a0154adeSJed Brown       // dim*dim contractions for grad
114a0154adeSJed Brown       for (CeedInt dim_1 = 0; dim_1 < BASIS_DIM; dim_1++) {
115a0154adeSJed Brown         CeedInt pre = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
116a0154adeSJed Brown         CeedInt post = 1;
117a0154adeSJed Brown         const CeedScalar *cur_u = u + elem * u_stride + dim_1 * u_dim_stride +
118a0154adeSJed Brown                                   comp * u_comp_stride;
119a0154adeSJed Brown         CeedScalar *cur_v = v + elem * v_stride + dim_1 * v_dim_stride + comp *
120a0154adeSJed Brown                             v_comp_stride;
121a0154adeSJed Brown         for (CeedInt dim_2 = 0; dim_2 < BASIS_DIM; dim_2++) {
122a0154adeSJed Brown           __syncthreads();
123a0154adeSJed Brown           // Update bufferfers used
124a0154adeSJed Brown           pre /= P;
125a0154adeSJed Brown           const CeedScalar *op = dim_1 == dim_2 ? s_grad_1d : s_interp_1d;
126a0154adeSJed Brown           const CeedScalar *in = dim_2 == 0
127a0154adeSJed Brown                                  ? cur_u
128a0154adeSJed Brown                                  : (dim_2 % 2 ? s_buffer_2 : s_buffer_1);
129a0154adeSJed Brown           CeedScalar *out = dim_2 == BASIS_DIM - 1
130a0154adeSJed Brown                             ? cur_v
131a0154adeSJed Brown                             : (dim_2 % 2 ? s_buffer_1 : s_buffer_2);
132a0154adeSJed Brown 
133a0154adeSJed Brown           // Contract along middle index
134a0154adeSJed Brown           const CeedInt writeLen = pre * post * Q;
135a0154adeSJed Brown           for (CeedInt k = i; k < writeLen; k += blockDim.x) {
136a0154adeSJed Brown             const CeedInt c = k % post;
137a0154adeSJed Brown             const CeedInt j = (k / post) % Q;
138a0154adeSJed Brown             const CeedInt a = k / (post * Q);
139a0154adeSJed Brown             CeedScalar v_k = 0;
140a0154adeSJed Brown             for (CeedInt b = 0; b < P; b++)
141a0154adeSJed Brown               v_k += op[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c];
142a0154adeSJed Brown 
143a0154adeSJed Brown             if (transpose && dim_2 == BASIS_DIM - 1)
144a0154adeSJed Brown               out[k] += v_k;
145a0154adeSJed Brown             else
146a0154adeSJed Brown               out[k] = v_k;
147a0154adeSJed Brown           }
148a0154adeSJed Brown 
149a0154adeSJed Brown           post *= Q;
150a0154adeSJed Brown         }
151a0154adeSJed Brown       }
152a0154adeSJed Brown     }
153a0154adeSJed Brown   }
154a0154adeSJed Brown }
155a0154adeSJed Brown 
156a0154adeSJed Brown //------------------------------------------------------------------------------
157a0154adeSJed Brown // 1D quadrature weights
158a0154adeSJed Brown //------------------------------------------------------------------------------
159a0154adeSJed Brown __device__ void Weight1d(const CeedInt num_elem, const CeedScalar *q_weight_1d,
160a0154adeSJed Brown                          CeedScalar *w) {
161a0154adeSJed Brown   const CeedInt i = threadIdx.x;
162a0154adeSJed Brown   if (i < BASIS_Q_1D) {
163a0154adeSJed Brown     const size_t elem = blockIdx.x;
164a0154adeSJed Brown     if (elem < num_elem)
165a0154adeSJed Brown       w[elem*BASIS_Q_1D + i] = q_weight_1d[i];
166a0154adeSJed Brown   }
167a0154adeSJed Brown }
168a0154adeSJed Brown 
169a0154adeSJed Brown //------------------------------------------------------------------------------
170a0154adeSJed Brown // 2D quadrature weights
171a0154adeSJed Brown //------------------------------------------------------------------------------
172a0154adeSJed Brown __device__ void Weight2d(const CeedInt num_elem, const CeedScalar *q_weight_1d,
173a0154adeSJed Brown                          CeedScalar *w) {
174a0154adeSJed Brown 
175a0154adeSJed Brown   const CeedInt i = threadIdx.x;
176a0154adeSJed Brown   const CeedInt j = threadIdx.y;
177a0154adeSJed Brown   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
178a0154adeSJed Brown     const size_t elem = blockIdx.x;
179a0154adeSJed Brown     if (elem < num_elem) {
180a0154adeSJed Brown       const size_t ind = (elem * BASIS_Q_1D + j) * BASIS_Q_1D + i;
181a0154adeSJed Brown       w[ind] = q_weight_1d[i] * q_weight_1d[j];
182a0154adeSJed Brown     }
183a0154adeSJed Brown   }
184a0154adeSJed Brown }
185a0154adeSJed Brown 
186a0154adeSJed Brown //------------------------------------------------------------------------------
187a0154adeSJed Brown // 3D quadrature weights
188a0154adeSJed Brown //------------------------------------------------------------------------------
189a0154adeSJed Brown __device__ void Weight3d(const CeedInt num_elem, const CeedScalar *q_weight_1d,
190a0154adeSJed Brown                          CeedScalar *w) {
191a0154adeSJed Brown   const CeedInt i = threadIdx.x;
192a0154adeSJed Brown   const CeedInt j = threadIdx.y;
193a0154adeSJed Brown   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
194a0154adeSJed Brown     const size_t elem = blockIdx.x;
195a0154adeSJed Brown     if (elem < num_elem) {
196a0154adeSJed Brown       for (CeedInt k = 0; k < BASIS_Q_1D; k++) {
197a0154adeSJed Brown         const size_t ind = ((elem * BASIS_Q_1D + k) * BASIS_Q_1D + j) * BASIS_Q_1D + i;
198a0154adeSJed Brown         w[ind] = q_weight_1d[i] * q_weight_1d[j] * q_weight_1d[k];
199a0154adeSJed Brown       }
200a0154adeSJed Brown     }
201a0154adeSJed Brown   }
202a0154adeSJed Brown }
203a0154adeSJed Brown 
204a0154adeSJed Brown //------------------------------------------------------------------------------
205a0154adeSJed Brown // Quadrature weights
206a0154adeSJed Brown //------------------------------------------------------------------------------
207a0154adeSJed Brown extern "C" __global__ void Weight(const CeedInt num_elem,
208a0154adeSJed Brown                                   const CeedScalar *__restrict__ q_weight_1d,
209a0154adeSJed Brown                                   CeedScalar *__restrict__ v) {
210a0154adeSJed Brown   if (BASIS_DIM == 1)
211a0154adeSJed Brown     Weight1d(num_elem, q_weight_1d, v);
212a0154adeSJed Brown   else if (BASIS_DIM == 2)
213a0154adeSJed Brown     Weight2d(num_elem, q_weight_1d, v);
214a0154adeSJed Brown   else if (BASIS_DIM == 3)
215a0154adeSJed Brown     Weight3d(num_elem, q_weight_1d, v);
216a0154adeSJed Brown }
217a0154adeSJed Brown 
218a0154adeSJed Brown //------------------------------------------------------------------------------
219