xref: /libCEED/include/ceed/jit-source/cuda/cuda-ref-basis-tensor.h (revision ca7355308a14805110a7d98dabf1f658d5ec16d5)
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 
8b2165e7aSSebastian Grimberg /// @file
9b2165e7aSSebastian Grimberg /// Internal header for CUDA tensor product basis
1094b7b29bSJeremy L Thompson #ifndef CEED_CUDA_REF_BASIS_TENSOR_H
1194b7b29bSJeremy L Thompson #define CEED_CUDA_REF_BASIS_TENSOR_H
12b2165e7aSSebastian Grimberg 
13c9c2c079SJeremy L Thompson #include <ceed.h>
14a0154adeSJed Brown 
15a0154adeSJed Brown //------------------------------------------------------------------------------
16a0154adeSJed Brown // Tensor Basis Kernels
17a0154adeSJed Brown //------------------------------------------------------------------------------
18a0154adeSJed Brown 
19a0154adeSJed Brown //------------------------------------------------------------------------------
20a0154adeSJed Brown // Interp
21a0154adeSJed Brown //------------------------------------------------------------------------------
222b730f8bSJeremy L Thompson extern "C" __global__ void Interp(const CeedInt num_elem, const CeedInt transpose, const CeedScalar *__restrict__ interp_1d,
232b730f8bSJeremy L Thompson                                   const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
24a0154adeSJed Brown   const CeedInt i = threadIdx.x;
25a0154adeSJed Brown 
26a0154adeSJed Brown   __shared__ CeedScalar s_mem[BASIS_Q_1D * BASIS_P_1D + 2 * BASIS_BUF_LEN];
27a0154adeSJed Brown   CeedScalar           *s_interp_1d = s_mem;
28a0154adeSJed Brown   CeedScalar           *s_buffer_1  = s_mem + BASIS_Q_1D * BASIS_P_1D;
29a0154adeSJed Brown   CeedScalar           *s_buffer_2  = s_buffer_1 + BASIS_BUF_LEN;
30a0154adeSJed Brown   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
31a0154adeSJed Brown     s_interp_1d[k] = interp_1d[k];
32a0154adeSJed Brown   }
33a0154adeSJed Brown 
34a0154adeSJed Brown   const CeedInt P             = transpose ? BASIS_Q_1D : BASIS_P_1D;
35a0154adeSJed Brown   const CeedInt Q             = transpose ? BASIS_P_1D : BASIS_Q_1D;
36a0154adeSJed Brown   const CeedInt stride_0      = transpose ? 1 : BASIS_P_1D;
37a0154adeSJed Brown   const CeedInt stride_1      = transpose ? BASIS_P_1D : 1;
38a0154adeSJed Brown   const CeedInt u_stride      = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
39a0154adeSJed Brown   const CeedInt v_stride      = transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
40a0154adeSJed Brown   const CeedInt u_comp_stride = num_elem * (transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
41a0154adeSJed Brown   const CeedInt v_comp_stride = num_elem * (transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
42a0154adeSJed Brown   const CeedInt u_size        = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
43a0154adeSJed Brown 
44a0154adeSJed Brown   // Apply basis element by element
45a0154adeSJed Brown   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
46a0154adeSJed Brown     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
47a0154adeSJed Brown       const CeedScalar *cur_u = u + elem * u_stride + comp * u_comp_stride;
48a0154adeSJed Brown       CeedScalar       *cur_v = v + elem * v_stride + comp * v_comp_stride;
49*ca735530SJeremy 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       CeedInt pre  = u_size;
54a0154adeSJed Brown       CeedInt post = 1;
55*ca735530SJeremy L Thompson 
56a0154adeSJed Brown       for (CeedInt d = 0; d < BASIS_DIM; d++) {
57a0154adeSJed Brown         __syncthreads();
582a86cc9dSSebastian Grimberg         // Update buffers used
59a0154adeSJed Brown         pre /= P;
60a0154adeSJed Brown         const CeedScalar *in  = d % 2 ? s_buffer_2 : s_buffer_1;
61a0154adeSJed Brown         CeedScalar       *out = d == BASIS_DIM - 1 ? cur_v : (d % 2 ? s_buffer_1 : s_buffer_2);
62a0154adeSJed Brown 
63a0154adeSJed Brown         // Contract along middle index
64a0154adeSJed Brown         const CeedInt writeLen = pre * post * Q;
65a0154adeSJed Brown         for (CeedInt k = i; k < writeLen; k += blockDim.x) {
66a0154adeSJed Brown           const CeedInt c = k % post;
67a0154adeSJed Brown           const CeedInt j = (k / post) % Q;
68a0154adeSJed Brown           const CeedInt a = k / (post * Q);
69a0154adeSJed Brown 
70a0154adeSJed Brown           CeedScalar vk = 0;
712b730f8bSJeremy 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];
72a0154adeSJed Brown           out[k] = vk;
73a0154adeSJed Brown         }
74a0154adeSJed Brown         post *= Q;
75a0154adeSJed Brown       }
76a0154adeSJed Brown     }
77a0154adeSJed Brown   }
78a0154adeSJed Brown }
79a0154adeSJed Brown 
80a0154adeSJed Brown //------------------------------------------------------------------------------
81a0154adeSJed Brown // Grad
82a0154adeSJed Brown //------------------------------------------------------------------------------
832b730f8bSJeremy L Thompson extern "C" __global__ void Grad(const CeedInt num_elem, const CeedInt transpose, const CeedScalar *__restrict__ interp_1d,
842b730f8bSJeremy L Thompson                                 const CeedScalar *__restrict__ grad_1d, const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) {
85a0154adeSJed Brown   const CeedInt i = threadIdx.x;
86a0154adeSJed Brown 
87a0154adeSJed Brown   __shared__ CeedScalar s_mem[2 * (BASIS_Q_1D * BASIS_P_1D + BASIS_BUF_LEN)];
88a0154adeSJed Brown   CeedScalar           *s_interp_1d = s_mem;
89a0154adeSJed Brown   CeedScalar           *s_grad_1d   = s_interp_1d + BASIS_Q_1D * BASIS_P_1D;
90a0154adeSJed Brown   CeedScalar           *s_buffer_1  = s_grad_1d + BASIS_Q_1D * BASIS_P_1D;
91a0154adeSJed Brown   CeedScalar           *s_buffer_2  = s_buffer_1 + BASIS_BUF_LEN;
92a0154adeSJed Brown   for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) {
93a0154adeSJed Brown     s_interp_1d[k] = interp_1d[k];
94a0154adeSJed Brown     s_grad_1d[k]   = grad_1d[k];
95a0154adeSJed Brown   }
96a0154adeSJed Brown 
97a0154adeSJed Brown   const CeedInt P             = transpose ? BASIS_Q_1D : BASIS_P_1D;
98a0154adeSJed Brown   const CeedInt Q             = transpose ? BASIS_P_1D : BASIS_Q_1D;
99a0154adeSJed Brown   const CeedInt stride_0      = transpose ? 1 : BASIS_P_1D;
100a0154adeSJed Brown   const CeedInt stride_1      = transpose ? BASIS_P_1D : 1;
101a0154adeSJed Brown   const CeedInt u_stride      = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
102a0154adeSJed Brown   const CeedInt v_stride      = transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS;
103a0154adeSJed Brown   const CeedInt u_comp_stride = num_elem * (transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES);
104a0154adeSJed Brown   const CeedInt v_comp_stride = num_elem * (transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS);
105a0154adeSJed Brown   const CeedInt u_dim_stride  = transpose ? num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP : 0;
106a0154adeSJed Brown   const CeedInt v_dim_stride  = transpose ? 0 : num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP;
107a0154adeSJed Brown 
108a0154adeSJed Brown   // Apply basis element by element
109a0154adeSJed Brown   for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) {
110a0154adeSJed Brown     for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) {
111a0154adeSJed Brown       // dim*dim contractions for grad
112a0154adeSJed Brown       for (CeedInt dim_1 = 0; dim_1 < BASIS_DIM; dim_1++) {
113a0154adeSJed Brown         CeedInt           pre   = transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES;
114a0154adeSJed Brown         CeedInt           post  = 1;
1152b730f8bSJeremy L Thompson         const CeedScalar *cur_u = u + elem * u_stride + dim_1 * u_dim_stride + comp * u_comp_stride;
1162b730f8bSJeremy L Thompson         CeedScalar       *cur_v = v + elem * v_stride + dim_1 * v_dim_stride + comp * v_comp_stride;
117*ca735530SJeremy L Thompson 
118a0154adeSJed Brown         for (CeedInt dim_2 = 0; dim_2 < BASIS_DIM; dim_2++) {
119a0154adeSJed Brown           __syncthreads();
120b2165e7aSSebastian Grimberg           // Update buffers used
121a0154adeSJed Brown           pre /= P;
122a0154adeSJed Brown           const CeedScalar *op  = dim_1 == dim_2 ? s_grad_1d : s_interp_1d;
1232b730f8bSJeremy L Thompson           const CeedScalar *in  = dim_2 == 0 ? cur_u : (dim_2 % 2 ? s_buffer_2 : s_buffer_1);
1242b730f8bSJeremy L Thompson           CeedScalar       *out = dim_2 == BASIS_DIM - 1 ? cur_v : (dim_2 % 2 ? s_buffer_1 : s_buffer_2);
125a0154adeSJed Brown 
126a0154adeSJed Brown           // Contract along middle index
127a0154adeSJed Brown           const CeedInt writeLen = pre * post * Q;
128a0154adeSJed Brown           for (CeedInt k = i; k < writeLen; k += blockDim.x) {
129a0154adeSJed Brown             const CeedInt c   = k % post;
130a0154adeSJed Brown             const CeedInt j   = (k / post) % Q;
131a0154adeSJed Brown             const CeedInt a   = k / (post * Q);
132a0154adeSJed Brown             CeedScalar    v_k = 0;
133*ca735530SJeremy L Thompson 
1342b730f8bSJeremy L Thompson             for (CeedInt b = 0; b < P; b++) v_k += op[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c];
135a0154adeSJed Brown 
1362b730f8bSJeremy L Thompson             if (transpose && dim_2 == BASIS_DIM - 1) out[k] += v_k;
1372b730f8bSJeremy L Thompson             else out[k] = v_k;
138a0154adeSJed Brown           }
139a0154adeSJed Brown           post *= Q;
140a0154adeSJed Brown         }
141a0154adeSJed Brown       }
142a0154adeSJed Brown     }
143a0154adeSJed Brown   }
144a0154adeSJed Brown }
145a0154adeSJed Brown 
146a0154adeSJed Brown //------------------------------------------------------------------------------
147a0154adeSJed Brown // 1D quadrature weights
148a0154adeSJed Brown //------------------------------------------------------------------------------
1492b730f8bSJeremy L Thompson __device__ void Weight1d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
150a0154adeSJed Brown   const CeedInt i = threadIdx.x;
151*ca735530SJeremy L Thompson 
152a0154adeSJed Brown   if (i < BASIS_Q_1D) {
153a0154adeSJed Brown     const size_t elem = blockIdx.x;
154*ca735530SJeremy L Thompson 
1552b730f8bSJeremy L Thompson     if (elem < num_elem) w[elem * BASIS_Q_1D + i] = q_weight_1d[i];
156a0154adeSJed Brown   }
157a0154adeSJed Brown }
158a0154adeSJed Brown 
159a0154adeSJed Brown //------------------------------------------------------------------------------
160a0154adeSJed Brown // 2D quadrature weights
161a0154adeSJed Brown //------------------------------------------------------------------------------
1622b730f8bSJeremy L Thompson __device__ void Weight2d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
163a0154adeSJed Brown   const CeedInt i = threadIdx.x;
164a0154adeSJed Brown   const CeedInt j = threadIdx.y;
165*ca735530SJeremy L Thompson 
166a0154adeSJed Brown   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
167a0154adeSJed Brown     const size_t elem = blockIdx.x;
168*ca735530SJeremy L Thompson 
169a0154adeSJed Brown     if (elem < num_elem) {
170a0154adeSJed Brown       const size_t ind = (elem * BASIS_Q_1D + j) * BASIS_Q_1D + i;
171*ca735530SJeremy L Thompson 
172a0154adeSJed Brown       w[ind] = q_weight_1d[i] * q_weight_1d[j];
173a0154adeSJed Brown     }
174a0154adeSJed Brown   }
175a0154adeSJed Brown }
176a0154adeSJed Brown 
177a0154adeSJed Brown //------------------------------------------------------------------------------
178a0154adeSJed Brown // 3D quadrature weights
179a0154adeSJed Brown //------------------------------------------------------------------------------
1802b730f8bSJeremy L Thompson __device__ void Weight3d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) {
181a0154adeSJed Brown   const CeedInt i = threadIdx.x;
182a0154adeSJed Brown   const CeedInt j = threadIdx.y;
183*ca735530SJeremy L Thompson 
184a0154adeSJed Brown   if (i < BASIS_Q_1D && j < BASIS_Q_1D) {
185a0154adeSJed Brown     const size_t elem = blockIdx.x;
186*ca735530SJeremy L Thompson 
187a0154adeSJed Brown     if (elem < num_elem) {
188a0154adeSJed Brown       for (CeedInt k = 0; k < BASIS_Q_1D; k++) {
189a0154adeSJed Brown         const size_t ind = ((elem * BASIS_Q_1D + k) * BASIS_Q_1D + j) * BASIS_Q_1D + i;
190*ca735530SJeremy L Thompson 
191a0154adeSJed Brown         w[ind] = q_weight_1d[i] * q_weight_1d[j] * q_weight_1d[k];
192a0154adeSJed Brown       }
193a0154adeSJed Brown     }
194a0154adeSJed Brown   }
195a0154adeSJed Brown }
196a0154adeSJed Brown 
197a0154adeSJed Brown //------------------------------------------------------------------------------
198a0154adeSJed Brown // Quadrature weights
199a0154adeSJed Brown //------------------------------------------------------------------------------
2002b730f8bSJeremy L Thompson extern "C" __global__ void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ v) {
2012b730f8bSJeremy L Thompson   if (BASIS_DIM == 1) Weight1d(num_elem, q_weight_1d, v);
2022b730f8bSJeremy L Thompson   else if (BASIS_DIM == 2) Weight2d(num_elem, q_weight_1d, v);
2032b730f8bSJeremy L Thompson   else if (BASIS_DIM == 3) Weight3d(num_elem, q_weight_1d, v);
204a0154adeSJed Brown }
205a0154adeSJed Brown 
206a0154adeSJed Brown //------------------------------------------------------------------------------
207b2165e7aSSebastian Grimberg 
20894b7b29bSJeremy L Thompson #endif  // CEED_CUDA_REF_BASIS_TENSOR_H
209