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