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 backend macro and type definitions for JiT source 10 #include <ceed/types.h> 11 12 //------------------------------------------------------------------------------ 13 // Load matrices for basis actions 14 //------------------------------------------------------------------------------ 15 template <int P, int Q> 16 inline __device__ void loadMatrix(SharedData_Cuda &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) { 17 for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i]; 18 } 19 20 //------------------------------------------------------------------------------ 21 // 1D 22 //------------------------------------------------------------------------------ 23 24 //------------------------------------------------------------------------------ 25 // L-vector -> E-vector, offsets provided 26 //------------------------------------------------------------------------------ 27 template <int NUM_COMP, int COMP_STRIDE, int P_1d> 28 inline __device__ void readDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 29 const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 30 if (data.t_id_x < P_1d) { 31 const CeedInt node = data.t_id_x; 32 const CeedInt ind = indices[node + elem * P_1d]; 33 34 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 35 } 36 } 37 38 //------------------------------------------------------------------------------ 39 // L-vector -> E-vector, strided 40 //------------------------------------------------------------------------------ 41 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 42 inline __device__ void readDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, 43 CeedScalar *__restrict__ r_u) { 44 if (data.t_id_x < P_1d) { 45 const CeedInt node = data.t_id_x; 46 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 47 48 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 49 } 50 } 51 52 //------------------------------------------------------------------------------ 53 // E-vector -> L-vector, offsets provided 54 //------------------------------------------------------------------------------ 55 template <int NUM_COMP, int COMP_STRIDE, int P_1d> 56 inline __device__ void writeDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 57 const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 58 if (data.t_id_x < P_1d) { 59 const CeedInt node = data.t_id_x; 60 const CeedInt ind = indices[node + elem * P_1d]; 61 62 for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 63 } 64 } 65 66 //------------------------------------------------------------------------------ 67 // E-vector -> L-vector, strided 68 //------------------------------------------------------------------------------ 69 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 70 inline __device__ void writeDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 71 CeedScalar *__restrict__ d_v) { 72 if (data.t_id_x < P_1d) { 73 const CeedInt node = data.t_id_x; 74 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 75 76 for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 77 } 78 } 79 80 //------------------------------------------------------------------------------ 81 // 2D 82 //------------------------------------------------------------------------------ 83 84 //------------------------------------------------------------------------------ 85 // L-vector -> E-vector, offsets provided 86 //------------------------------------------------------------------------------ 87 template <int NUM_COMP, int COMP_STRIDE, int P_1d> 88 inline __device__ void readDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 89 const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 90 if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 91 const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 92 const CeedInt ind = indices[node + elem * P_1d * P_1d]; 93 94 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 95 } 96 } 97 98 //------------------------------------------------------------------------------ 99 // L-vector -> E-vector, strided 100 //------------------------------------------------------------------------------ 101 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 102 inline __device__ void readDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, 103 CeedScalar *__restrict__ r_u) { 104 if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 105 const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 106 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 107 108 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 109 } 110 } 111 112 //------------------------------------------------------------------------------ 113 // E-vector -> L-vector, offsets provided 114 //------------------------------------------------------------------------------ 115 template <int NUM_COMP, int COMP_STRIDE, int P_1d> 116 inline __device__ void writeDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 117 const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 118 if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 119 const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 120 const CeedInt ind = indices[node + elem * P_1d * P_1d]; 121 122 for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 123 } 124 } 125 126 //------------------------------------------------------------------------------ 127 // E-vector -> L-vector, strided 128 //------------------------------------------------------------------------------ 129 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 130 inline __device__ void writeDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 131 CeedScalar *__restrict__ d_v) { 132 if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 133 const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 134 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 135 136 for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 137 } 138 } 139 140 //------------------------------------------------------------------------------ 141 // 3D 142 //------------------------------------------------------------------------------ 143 144 //------------------------------------------------------------------------------ 145 // L-vector -> E-vector, offsets provided 146 //------------------------------------------------------------------------------ 147 // TODO: remove "Dofs" and "Quads" in the following function names? 148 // - readDofsOffset3d -> readOffset3d ? 149 // - readDofsStrided3d -> readStrided3d ? 150 // - readSliceQuadsOffset3d -> readSliceOffset3d ? 151 // - readSliceQuadsStrided3d -> readSliceStrided3d ? 152 // - writeDofsOffset3d -> writeOffset3d ? 153 // - writeDofsStrided3d -> writeStrided3d ? 154 template <int NUM_COMP, int COMP_STRIDE, int P_1d> 155 inline __device__ void readDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 156 const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 157 if (data.t_id_x < P_1d && data.t_id_y < P_1d) 158 for (CeedInt z = 0; z < P_1d; z++) { 159 const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 160 const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 161 162 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + COMP_STRIDE * comp]; 163 } 164 } 165 166 //------------------------------------------------------------------------------ 167 // L-vector -> E-vector, strided 168 //------------------------------------------------------------------------------ 169 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 170 inline __device__ void readDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, 171 CeedScalar *__restrict__ r_u) { 172 if (data.t_id_x < P_1d && data.t_id_y < P_1d) 173 for (CeedInt z = 0; z < P_1d; z++) { 174 const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 175 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 176 177 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP]; 178 } 179 } 180 181 //------------------------------------------------------------------------------ 182 // E-vector -> Q-vector, offests provided 183 //------------------------------------------------------------------------------ 184 template <int NUM_COMP, int COMP_STRIDE, int Q_1d> 185 inline __device__ void readSliceQuadsOffset3d(SharedData_Cuda &data, const CeedInt nquads, const CeedInt elem, const CeedInt q, 186 const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 187 if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 188 const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 189 const CeedInt ind = indices[node + elem * Q_1d * Q_1d * Q_1d]; 190 191 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 192 } 193 } 194 195 //------------------------------------------------------------------------------ 196 // E-vector -> Q-vector, strided 197 //------------------------------------------------------------------------------ 198 template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 199 inline __device__ void readSliceQuadsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u, 200 CeedScalar *__restrict__ r_u) { 201 if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 202 const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 203 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 204 205 for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 206 } 207 } 208 209 //------------------------------------------------------------------------------ 210 // E-vector -> L-vector, offsets provided 211 //------------------------------------------------------------------------------ 212 template <int NUM_COMP, int COMP_STRIDE, int P_1d> 213 inline __device__ void writeDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 214 const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 215 if (data.t_id_x < P_1d && data.t_id_y < P_1d) 216 for (CeedInt z = 0; z < P_1d; z++) { 217 const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 218 const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 219 220 for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]); 221 } 222 } 223 224 //------------------------------------------------------------------------------ 225 // E-vector -> L-vector, strided 226 //------------------------------------------------------------------------------ 227 template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 228 inline __device__ void writeDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 229 CeedScalar *__restrict__ d_v) { 230 if (data.t_id_x < P_1d && data.t_id_y < P_1d) 231 for (CeedInt z = 0; z < P_1d; z++) { 232 const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 233 const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 234 235 for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d]; 236 } 237 } 238 239 //------------------------------------------------------------------------------ 240 // 3D collocated derivatives computation 241 //------------------------------------------------------------------------------ 242 template <int NUM_COMP, int Q_1d> 243 inline __device__ void gradCollo3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 244 CeedScalar *__restrict__ r_V) { 245 if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 246 for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 247 data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d]; 248 __syncthreads(); 249 // X derivative 250 r_V[comp + 0 * NUM_COMP] = 0.0; 251 for (CeedInt i = 0; i < Q_1d; i++) 252 r_V[comp + 0 * NUM_COMP] += c_G[i + data.t_id_x * Q_1d] * data.slice[i + data.t_id_y * T_1D]; // Contract x direction (X derivative) 253 // Y derivative 254 r_V[comp + 1 * NUM_COMP] = 0.0; 255 for (CeedInt i = 0; i < Q_1d; i++) 256 r_V[comp + 1 * NUM_COMP] += c_G[i + data.t_id_y * Q_1d] * data.slice[data.t_id_x + i * T_1D]; // Contract y direction (Y derivative) 257 // Z derivative 258 r_V[comp + 2 * NUM_COMP] = 0.0; 259 for (CeedInt i = 0; i < Q_1d; i++) r_V[comp + 2 * NUM_COMP] += c_G[i + q * Q_1d] * r_U[i + comp * Q_1d]; // Contract z direction (Z derivative) 260 __syncthreads(); 261 } 262 } 263 } 264 265 //------------------------------------------------------------------------------ 266 // 3D collocated derivatives transpose 267 //------------------------------------------------------------------------------ 268 template <int NUM_COMP, int Q_1d> 269 inline __device__ void gradColloTranspose3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 270 CeedScalar *__restrict__ r_V) { 271 if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 272 for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 273 // X derivative 274 data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP]; 275 __syncthreads(); 276 for (CeedInt i = 0; i < Q_1d; i++) 277 r_V[q + comp * Q_1d] += c_G[data.t_id_x + i * Q_1d] * data.slice[i + data.t_id_y * T_1D]; // Contract x direction (X derivative) 278 __syncthreads(); 279 // Y derivative 280 data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP]; 281 __syncthreads(); 282 for (CeedInt i = 0; i < Q_1d; i++) 283 r_V[q + comp * Q_1d] += c_G[data.t_id_y + i * Q_1d] * data.slice[data.t_id_x + i * T_1D]; // Contract y direction (Y derivative) 284 __syncthreads(); 285 // Z derivative 286 for (CeedInt i = 0; i < Q_1d; i++) 287 r_V[i + comp * Q_1d] += c_G[i + q * Q_1d] * r_U[comp + 2 * NUM_COMP]; // PARTIAL contract z direction (Z derivative) 288 } 289 } 290 } 291