1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 10c0b5abf0SJeremy L Thompson #include <ceed/types.h> 11a0154adeSJed Brown 12a0154adeSJed Brown //------------------------------------------------------------------------------ 13a0154adeSJed Brown // Tensor Basis Kernels 14a0154adeSJed Brown //------------------------------------------------------------------------------ 15a0154adeSJed Brown 16a0154adeSJed Brown //------------------------------------------------------------------------------ 17a0154adeSJed Brown // Interp 18a0154adeSJed Brown //------------------------------------------------------------------------------ 190b2e4913SJeremy L Thompson extern "C" __global__ void Interp(const CeedInt num_elem, const CeedInt is_transpose, const CeedScalar *__restrict__ interp_1d, 202b730f8bSJeremy L Thompson const CeedScalar *__restrict__ u, 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 310b2e4913SJeremy L Thompson const CeedInt P = is_transpose ? BASIS_Q_1D : BASIS_P_1D; 320b2e4913SJeremy L Thompson const CeedInt Q = is_transpose ? BASIS_P_1D : BASIS_Q_1D; 330b2e4913SJeremy L Thompson const CeedInt stride_0 = is_transpose ? 1 : BASIS_P_1D; 340b2e4913SJeremy L Thompson const CeedInt stride_1 = is_transpose ? BASIS_P_1D : 1; 350b2e4913SJeremy L Thompson const CeedInt u_stride = is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES; 360b2e4913SJeremy L Thompson const CeedInt v_stride = is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS; 370b2e4913SJeremy L Thompson const CeedInt u_comp_stride = num_elem * (is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES); 380b2e4913SJeremy L Thompson const CeedInt v_comp_stride = num_elem * (is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS); 390b2e4913SJeremy L Thompson const CeedInt u_size = is_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++) { 44db2becc9SJeremy L Thompson const CeedScalar *cur_u = &u[elem * u_stride + comp * u_comp_stride]; 45db2becc9SJeremy L Thompson CeedScalar *cur_v = &v[elem * v_stride + comp * v_comp_stride]; 46672b0f2aSSebastian Grimberg CeedInt pre = u_size; 47672b0f2aSSebastian Grimberg CeedInt post = 1; 48ca735530SJeremy L Thompson 49a0154adeSJed Brown for (CeedInt d = 0; d < BASIS_DIM; d++) { 50a0154adeSJed Brown __syncthreads(); 512a86cc9dSSebastian Grimberg // Update buffers used 52a0154adeSJed Brown pre /= P; 532d903c70SJeremy L Thompson const CeedScalar *in = d == 0 ? cur_u : (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); 55672b0f2aSSebastian Grimberg const CeedInt writeLen = pre * post * Q; 56a0154adeSJed Brown 57a0154adeSJed Brown // Contract along middle index 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); 62db2becc9SJeremy L Thompson CeedScalar v_k = 0; 63672b0f2aSSebastian Grimberg 64db2becc9SJeremy L Thompson for (CeedInt b = 0; b < P; b++) v_k += s_interp_1d[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c]; 65db2becc9SJeremy L Thompson if (is_transpose && d == BASIS_DIM - 1) out[k] += v_k; 66db2becc9SJeremy L Thompson else out[k] = v_k; 67a0154adeSJed Brown } 68a0154adeSJed Brown post *= Q; 69a0154adeSJed Brown } 70a0154adeSJed Brown } 71a0154adeSJed Brown } 72a0154adeSJed Brown } 73a0154adeSJed Brown 74a0154adeSJed Brown //------------------------------------------------------------------------------ 75a0154adeSJed Brown // Grad 76a0154adeSJed Brown //------------------------------------------------------------------------------ 770b2e4913SJeremy L Thompson extern "C" __global__ void Grad(const CeedInt num_elem, const CeedInt is_transpose, const CeedScalar *__restrict__ interp_1d, 782b730f8bSJeremy L Thompson const CeedScalar *__restrict__ grad_1d, const CeedScalar *__restrict__ u, CeedScalar *__restrict__ v) { 79a0154adeSJed Brown const CeedInt i = threadIdx.x; 80a0154adeSJed Brown 81a0154adeSJed Brown __shared__ CeedScalar s_mem[2 * (BASIS_Q_1D * BASIS_P_1D + BASIS_BUF_LEN)]; 82a0154adeSJed Brown CeedScalar *s_interp_1d = s_mem; 83a0154adeSJed Brown CeedScalar *s_grad_1d = s_interp_1d + BASIS_Q_1D * BASIS_P_1D; 84a0154adeSJed Brown CeedScalar *s_buffer_1 = s_grad_1d + BASIS_Q_1D * BASIS_P_1D; 85a0154adeSJed Brown CeedScalar *s_buffer_2 = s_buffer_1 + BASIS_BUF_LEN; 86a0154adeSJed Brown for (CeedInt k = i; k < BASIS_Q_1D * BASIS_P_1D; k += blockDim.x) { 87a0154adeSJed Brown s_interp_1d[k] = interp_1d[k]; 88a0154adeSJed Brown s_grad_1d[k] = grad_1d[k]; 89a0154adeSJed Brown } 90a0154adeSJed Brown 910b2e4913SJeremy L Thompson const CeedInt P = is_transpose ? BASIS_Q_1D : BASIS_P_1D; 920b2e4913SJeremy L Thompson const CeedInt Q = is_transpose ? BASIS_P_1D : BASIS_Q_1D; 930b2e4913SJeremy L Thompson const CeedInt stride_0 = is_transpose ? 1 : BASIS_P_1D; 940b2e4913SJeremy L Thompson const CeedInt stride_1 = is_transpose ? BASIS_P_1D : 1; 950b2e4913SJeremy L Thompson const CeedInt u_stride = is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES; 960b2e4913SJeremy L Thompson const CeedInt v_stride = is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS; 970b2e4913SJeremy L Thompson const CeedInt u_comp_stride = num_elem * (is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES); 980b2e4913SJeremy L Thompson const CeedInt v_comp_stride = num_elem * (is_transpose ? BASIS_NUM_NODES : BASIS_NUM_QPTS); 990b2e4913SJeremy L Thompson const CeedInt u_dim_stride = is_transpose ? num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP : 0; 1000b2e4913SJeremy L Thompson const CeedInt v_dim_stride = is_transpose ? 0 : num_elem * BASIS_NUM_QPTS * BASIS_NUM_COMP; 101a0154adeSJed Brown 102a0154adeSJed Brown // Apply basis element by element 103a0154adeSJed Brown for (CeedInt elem = blockIdx.x; elem < num_elem; elem += gridDim.x) { 104a0154adeSJed Brown for (CeedInt comp = 0; comp < BASIS_NUM_COMP; comp++) { 105a0154adeSJed Brown // dim*dim contractions for grad 106a0154adeSJed Brown for (CeedInt dim_1 = 0; dim_1 < BASIS_DIM; dim_1++) { 1070b2e4913SJeremy L Thompson CeedInt pre = is_transpose ? BASIS_NUM_QPTS : BASIS_NUM_NODES; 108a0154adeSJed Brown CeedInt post = 1; 109db2becc9SJeremy L Thompson const CeedScalar *cur_u = &u[elem * u_stride + dim_1 * u_dim_stride + comp * u_comp_stride]; 110db2becc9SJeremy L Thompson CeedScalar *cur_v = &v[elem * v_stride + dim_1 * v_dim_stride + comp * v_comp_stride]; 111ca735530SJeremy L Thompson 112a0154adeSJed Brown for (CeedInt dim_2 = 0; dim_2 < BASIS_DIM; dim_2++) { 113a0154adeSJed Brown __syncthreads(); 114b2165e7aSSebastian Grimberg // Update buffers used 115a0154adeSJed Brown pre /= P; 116a0154adeSJed Brown const CeedScalar *op = dim_1 == dim_2 ? s_grad_1d : s_interp_1d; 1172b730f8bSJeremy L Thompson const CeedScalar *in = dim_2 == 0 ? cur_u : (dim_2 % 2 ? s_buffer_2 : s_buffer_1); 1182b730f8bSJeremy L Thompson CeedScalar *out = dim_2 == BASIS_DIM - 1 ? cur_v : (dim_2 % 2 ? s_buffer_1 : s_buffer_2); 119672b0f2aSSebastian Grimberg const CeedInt writeLen = pre * post * Q; 120a0154adeSJed Brown 121a0154adeSJed Brown // Contract along middle index 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; 127ca735530SJeremy L Thompson 1282b730f8bSJeremy L Thompson for (CeedInt b = 0; b < P; b++) v_k += op[j * stride_0 + b * stride_1] * in[(a * P + b) * post + c]; 1290b2e4913SJeremy L Thompson if (is_transpose && dim_2 == BASIS_DIM - 1) out[k] += v_k; 1302b730f8bSJeremy L Thompson else out[k] = v_k; 131a0154adeSJed Brown } 132a0154adeSJed Brown post *= Q; 133a0154adeSJed Brown } 134a0154adeSJed Brown } 135a0154adeSJed Brown } 136a0154adeSJed Brown } 137a0154adeSJed Brown } 138a0154adeSJed Brown 139a0154adeSJed Brown //------------------------------------------------------------------------------ 140a0154adeSJed Brown // 1D quadrature weights 141a0154adeSJed Brown //------------------------------------------------------------------------------ 1422b730f8bSJeremy L Thompson __device__ void Weight1d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) { 143a0154adeSJed Brown const CeedInt i = threadIdx.x; 144ca735530SJeremy L Thompson 145a0154adeSJed Brown if (i < BASIS_Q_1D) { 146a0154adeSJed Brown const size_t elem = blockIdx.x; 147ca735530SJeremy L Thompson 1482b730f8bSJeremy L Thompson if (elem < num_elem) w[elem * BASIS_Q_1D + i] = q_weight_1d[i]; 149a0154adeSJed Brown } 150a0154adeSJed Brown } 151a0154adeSJed Brown 152a0154adeSJed Brown //------------------------------------------------------------------------------ 153a0154adeSJed Brown // 2D quadrature weights 154a0154adeSJed Brown //------------------------------------------------------------------------------ 1552b730f8bSJeremy L Thompson __device__ void Weight2d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) { 156a0154adeSJed Brown const CeedInt i = threadIdx.x; 157a0154adeSJed Brown const CeedInt j = threadIdx.y; 158ca735530SJeremy L Thompson 159a0154adeSJed Brown if (i < BASIS_Q_1D && j < BASIS_Q_1D) { 160a0154adeSJed Brown const size_t elem = blockIdx.x; 161ca735530SJeremy L Thompson 162a0154adeSJed Brown if (elem < num_elem) { 163a0154adeSJed Brown const size_t ind = (elem * BASIS_Q_1D + j) * BASIS_Q_1D + i; 164ca735530SJeremy L Thompson 165a0154adeSJed Brown w[ind] = q_weight_1d[i] * q_weight_1d[j]; 166a0154adeSJed Brown } 167a0154adeSJed Brown } 168a0154adeSJed Brown } 169a0154adeSJed Brown 170a0154adeSJed Brown //------------------------------------------------------------------------------ 171a0154adeSJed Brown // 3D quadrature weights 172a0154adeSJed Brown //------------------------------------------------------------------------------ 1732b730f8bSJeremy L Thompson __device__ void Weight3d(const CeedInt num_elem, const CeedScalar *q_weight_1d, CeedScalar *w) { 174a0154adeSJed Brown const CeedInt i = threadIdx.x; 175a0154adeSJed Brown const CeedInt j = threadIdx.y; 176ca735530SJeremy L Thompson 177a0154adeSJed Brown if (i < BASIS_Q_1D && j < BASIS_Q_1D) { 178a0154adeSJed Brown const size_t elem = blockIdx.x; 179ca735530SJeremy L Thompson 180a0154adeSJed Brown if (elem < num_elem) { 181a0154adeSJed Brown for (CeedInt k = 0; k < BASIS_Q_1D; k++) { 182a0154adeSJed Brown const size_t ind = ((elem * BASIS_Q_1D + k) * BASIS_Q_1D + j) * BASIS_Q_1D + i; 183ca735530SJeremy L Thompson 184a0154adeSJed Brown w[ind] = q_weight_1d[i] * q_weight_1d[j] * q_weight_1d[k]; 185a0154adeSJed Brown } 186a0154adeSJed Brown } 187a0154adeSJed Brown } 188a0154adeSJed Brown } 189a0154adeSJed Brown 190a0154adeSJed Brown //------------------------------------------------------------------------------ 191a0154adeSJed Brown // Quadrature weights 192a0154adeSJed Brown //------------------------------------------------------------------------------ 1932b730f8bSJeremy L Thompson extern "C" __global__ void Weight(const CeedInt num_elem, const CeedScalar *__restrict__ q_weight_1d, CeedScalar *__restrict__ v) { 1942b730f8bSJeremy L Thompson if (BASIS_DIM == 1) Weight1d(num_elem, q_weight_1d, v); 1952b730f8bSJeremy L Thompson else if (BASIS_DIM == 2) Weight2d(num_elem, q_weight_1d, v); 1962b730f8bSJeremy L Thompson else if (BASIS_DIM == 3) Weight3d(num_elem, q_weight_1d, v); 197a0154adeSJed Brown } 198