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