1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 29e201c85SYohann // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 39e201c85SYohann // 49e201c85SYohann // SPDX-License-Identifier: BSD-2-Clause 59e201c85SYohann // 69e201c85SYohann // This file is part of CEED: http://github.com/ceed 79e201c85SYohann 89e201c85SYohann /// @file 99e201c85SYohann /// Internal header for HIP backend macro and type definitions for JiT source 1094b7b29bSJeremy L Thompson #ifndef CEED_HIP_GEN_TEMPLATES_H 1194b7b29bSJeremy L Thompson #define CEED_HIP_GEN_TEMPLATES_H 129e201c85SYohann 131b7492f8SSebastian Grimberg #include <ceed.h> 149e201c85SYohann 159e201c85SYohann //------------------------------------------------------------------------------ 169e201c85SYohann // Load matrices for basis actions 179e201c85SYohann //------------------------------------------------------------------------------ 189e201c85SYohann template <int P, int Q> 19672b0f2aSSebastian Grimberg inline __device__ void loadMatrix(SharedData_Hip &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) { 202b730f8bSJeremy L Thompson for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i]; 219e201c85SYohann } 229e201c85SYohann 239e201c85SYohann //------------------------------------------------------------------------------ 249e201c85SYohann // 1D 259e201c85SYohann //------------------------------------------------------------------------------ 269e201c85SYohann 279e201c85SYohann //------------------------------------------------------------------------------ 289e201c85SYohann // L-vector -> E-vector, offsets provided 299e201c85SYohann //------------------------------------------------------------------------------ 30672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 31672b0f2aSSebastian Grimberg inline __device__ void readDofsOffset1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 32672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 33672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 349e201c85SYohann const CeedInt node = data.t_id_x; 35672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d]; 36672b0f2aSSebastian Grimberg 37672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 389e201c85SYohann } 399e201c85SYohann } 409e201c85SYohann 419e201c85SYohann //------------------------------------------------------------------------------ 429e201c85SYohann // L-vector -> E-vector, strided 439e201c85SYohann //------------------------------------------------------------------------------ 44672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 45672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 46672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 479e201c85SYohann const CeedInt node = data.t_id_x; 489e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 49672b0f2aSSebastian Grimberg 50672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 519e201c85SYohann } 529e201c85SYohann } 539e201c85SYohann 549e201c85SYohann //------------------------------------------------------------------------------ 559e201c85SYohann // E-vector -> L-vector, offsets provided 569e201c85SYohann //------------------------------------------------------------------------------ 57672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 58672b0f2aSSebastian Grimberg inline __device__ void writeDofsOffset1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 59672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 60672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 619e201c85SYohann const CeedInt node = data.t_id_x; 62672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d]; 63672b0f2aSSebastian Grimberg 64672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 659e201c85SYohann } 669e201c85SYohann } 679e201c85SYohann 689e201c85SYohann //------------------------------------------------------------------------------ 699e201c85SYohann // E-vector -> L-vector, strided 709e201c85SYohann //------------------------------------------------------------------------------ 71672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 72672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 73672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 74672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 759e201c85SYohann const CeedInt node = data.t_id_x; 769e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 77672b0f2aSSebastian Grimberg 78672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 799e201c85SYohann } 809e201c85SYohann } 819e201c85SYohann 829e201c85SYohann //------------------------------------------------------------------------------ 839e201c85SYohann // 2D 849e201c85SYohann //------------------------------------------------------------------------------ 859e201c85SYohann 869e201c85SYohann //------------------------------------------------------------------------------ 879e201c85SYohann // L-vector -> E-vector, offsets provided 889e201c85SYohann //------------------------------------------------------------------------------ 89672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 90672b0f2aSSebastian Grimberg inline __device__ void readDofsOffset2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 91672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 92672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 93672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 94672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d]; 95672b0f2aSSebastian Grimberg 96672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 979e201c85SYohann } 989e201c85SYohann } 999e201c85SYohann 1009e201c85SYohann //------------------------------------------------------------------------------ 1019e201c85SYohann // L-vector -> E-vector, strided 1029e201c85SYohann //------------------------------------------------------------------------------ 103672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 104672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 105672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 106672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 1079e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 108672b0f2aSSebastian Grimberg 109672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 1109e201c85SYohann } 1119e201c85SYohann } 1129e201c85SYohann 1139e201c85SYohann //------------------------------------------------------------------------------ 1149e201c85SYohann // E-vector -> L-vector, offsets provided 1159e201c85SYohann //------------------------------------------------------------------------------ 116672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 117672b0f2aSSebastian Grimberg inline __device__ void writeDofsOffset2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 118672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 119672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 120672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 121672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d]; 122672b0f2aSSebastian Grimberg 123672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 1249e201c85SYohann } 1259e201c85SYohann } 1269e201c85SYohann 1279e201c85SYohann //------------------------------------------------------------------------------ 1289e201c85SYohann // E-vector -> L-vector, strided 1299e201c85SYohann //------------------------------------------------------------------------------ 130672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 131672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 132672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 133672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 134672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 1359e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 136672b0f2aSSebastian Grimberg 137672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 1389e201c85SYohann } 1399e201c85SYohann } 1409e201c85SYohann 1419e201c85SYohann //------------------------------------------------------------------------------ 1429e201c85SYohann // 3D 1439e201c85SYohann //------------------------------------------------------------------------------ 1449e201c85SYohann 1459e201c85SYohann //------------------------------------------------------------------------------ 1469e201c85SYohann // L-vector -> E-vector, offsets provided 1479e201c85SYohann //------------------------------------------------------------------------------ 148672b0f2aSSebastian Grimberg // TODO: remove "Dofs" and "Quads" in the following function names? 149672b0f2aSSebastian Grimberg // - readDofsOffset3d -> readOffset3d ? 150672b0f2aSSebastian Grimberg // - readDofsStrided3d -> readStrided3d ? 151672b0f2aSSebastian Grimberg // - readSliceQuadsOffset3d -> readSliceOffset3d ? 152672b0f2aSSebastian Grimberg // - readSliceQuadsStrided3d -> readSliceStrided3d ? 153672b0f2aSSebastian Grimberg // - writeDofsOffset3d -> writeOffset3d ? 154672b0f2aSSebastian Grimberg // - writeDofsStrided3d -> writeStrided3d ? 155672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 156672b0f2aSSebastian Grimberg inline __device__ void readDofsOffset3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 157672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 158672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 159672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 160672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 161672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 162672b0f2aSSebastian Grimberg 163672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + COMP_STRIDE * comp]; 1649e201c85SYohann } 1659e201c85SYohann } 1669e201c85SYohann 1679e201c85SYohann //------------------------------------------------------------------------------ 1689e201c85SYohann // L-vector -> E-vector, strided 1699e201c85SYohann //------------------------------------------------------------------------------ 170672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 171672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 172672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 173672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 174672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 1759e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 176672b0f2aSSebastian Grimberg 177672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP]; 1789e201c85SYohann } 1799e201c85SYohann } 1809e201c85SYohann 1819e201c85SYohann //------------------------------------------------------------------------------ 1829e201c85SYohann // E-vector -> Q-vector, offests provided 1839e201c85SYohann //------------------------------------------------------------------------------ 184672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int Q_1d> 185672b0f2aSSebastian Grimberg inline __device__ void readSliceQuadsOffset3d(SharedData_Hip &data, const CeedInt nquads, const CeedInt elem, const CeedInt q, 186672b0f2aSSebastian Grimberg const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 187672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 188672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 189672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * Q_1d * Q_1d * Q_1d]; 190672b0f2aSSebastian Grimberg 191672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 1929e201c85SYohann } 1939e201c85SYohann } 1949e201c85SYohann 1959e201c85SYohann //------------------------------------------------------------------------------ 1969e201c85SYohann // E-vector -> Q-vector, strided 1979e201c85SYohann //------------------------------------------------------------------------------ 198672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 199672b0f2aSSebastian Grimberg inline __device__ void readSliceQuadsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u, 200672b0f2aSSebastian Grimberg CeedScalar *__restrict__ r_u) { 201672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 202672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 2039e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 204672b0f2aSSebastian Grimberg 205672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 2069e201c85SYohann } 2079e201c85SYohann } 2089e201c85SYohann 2099e201c85SYohann //------------------------------------------------------------------------------ 2109e201c85SYohann // E-vector -> L-vector, offsets provided 2119e201c85SYohann //------------------------------------------------------------------------------ 212672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 213672b0f2aSSebastian Grimberg inline __device__ void writeDofsOffset3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 214672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 215672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 216672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 217672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 218672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 219672b0f2aSSebastian Grimberg 220672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]); 2219e201c85SYohann } 2229e201c85SYohann } 2239e201c85SYohann 2249e201c85SYohann //------------------------------------------------------------------------------ 2259e201c85SYohann // E-vector -> L-vector, strided 2269e201c85SYohann //------------------------------------------------------------------------------ 227672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 228672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 229672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 230672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 231672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 232672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 2339e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 234672b0f2aSSebastian Grimberg 235672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d]; 2369e201c85SYohann } 2379e201c85SYohann } 2389e201c85SYohann 2399e201c85SYohann //------------------------------------------------------------------------------ 2409e201c85SYohann // 3D collocated derivatives computation 2419e201c85SYohann //------------------------------------------------------------------------------ 242672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d> 2432b730f8bSJeremy L Thompson inline __device__ void gradCollo3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 2442b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 245672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 246672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 247672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d]; 2489e201c85SYohann __syncthreads(); 2499e201c85SYohann // X derivative 250672b0f2aSSebastian Grimberg r_V[comp + 0 * NUM_COMP] = 0.0; 251672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 252672b0f2aSSebastian Grimberg 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) 2539e201c85SYohann // Y derivative 254672b0f2aSSebastian Grimberg r_V[comp + 1 * NUM_COMP] = 0.0; 255672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 256672b0f2aSSebastian Grimberg 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) 2579e201c85SYohann // Z derivative 258672b0f2aSSebastian Grimberg r_V[comp + 2 * NUM_COMP] = 0.0; 259672b0f2aSSebastian Grimberg 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) 2609e201c85SYohann __syncthreads(); 2619e201c85SYohann } 2629e201c85SYohann } 2639e201c85SYohann } 2649e201c85SYohann 2659e201c85SYohann //------------------------------------------------------------------------------ 2669e201c85SYohann // 3D collocated derivatives transpose 2679e201c85SYohann //------------------------------------------------------------------------------ 268672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d> 2692b730f8bSJeremy L Thompson inline __device__ void gradColloTranspose3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 2702b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 271672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 272672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 2739e201c85SYohann // X derivative 274672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP]; 2759e201c85SYohann __syncthreads(); 276672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 277672b0f2aSSebastian Grimberg 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) 2789e201c85SYohann __syncthreads(); 2799e201c85SYohann // Y derivative 280672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP]; 2819e201c85SYohann __syncthreads(); 282672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 283672b0f2aSSebastian Grimberg 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) 2849e201c85SYohann __syncthreads(); 2859e201c85SYohann // Z derivative 286672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 287672b0f2aSSebastian Grimberg r_V[i + comp * Q_1d] += c_G[i + q * Q_1d] * r_U[comp + 2 * NUM_COMP]; // PARTIAL contract z direction (Z derivative) 2889e201c85SYohann } 2899e201c85SYohann } 2909e201c85SYohann } 2919e201c85SYohann 29294b7b29bSJeremy L Thompson #endif // CEED_HIP_GEN_TEMPLATES_H 293