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