1d275d636SJeremy L Thompson // Copyright (c) 2017-2025, 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 10c0b5abf0SJeremy L Thompson #include <ceed/types.h> 119e201c85SYohann 129e201c85SYohann //------------------------------------------------------------------------------ 139e201c85SYohann // Load matrices for basis actions 149e201c85SYohann //------------------------------------------------------------------------------ 159e201c85SYohann template <int P, int Q> 16f815fac9SJeremy L Thompson inline __device__ void LoadMatrix(SharedData_Hip &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) { 172b730f8bSJeremy L Thompson for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i]; 189e201c85SYohann } 199e201c85SYohann 209e201c85SYohann //------------------------------------------------------------------------------ 213a2968d6SJeremy L Thompson // AtPoints 223a2968d6SJeremy L Thompson //------------------------------------------------------------------------------ 233a2968d6SJeremy L Thompson 243a2968d6SJeremy L Thompson //------------------------------------------------------------------------------ 253a2968d6SJeremy L Thompson // L-vector -> single point 263a2968d6SJeremy L Thompson //------------------------------------------------------------------------------ 273a2968d6SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int NUM_PTS> 283a2968d6SJeremy L Thompson inline __device__ void ReadPoint(SharedData_Hip &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem, 293a2968d6SJeremy L Thompson const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) { 303a2968d6SJeremy L Thompson const CeedInt ind = indices[p + elem * NUM_PTS]; 313a2968d6SJeremy L Thompson 323a2968d6SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 333a2968d6SJeremy L Thompson r_u[comp] = d_u[ind + comp * COMP_STRIDE]; 343a2968d6SJeremy L Thompson } 353a2968d6SJeremy L Thompson } 363a2968d6SJeremy L Thompson 373a2968d6SJeremy L Thompson //------------------------------------------------------------------------------ 383a2968d6SJeremy L Thompson // Single point -> L-vector 393a2968d6SJeremy L Thompson //------------------------------------------------------------------------------ 403a2968d6SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int NUM_PTS> 413a2968d6SJeremy L Thompson inline __device__ void WritePoint(SharedData_Hip &data, const CeedInt elem, const CeedInt p, const CeedInt points_in_elem, 423a2968d6SJeremy L Thompson const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ r_u, CeedScalar *d_u) { 433a2968d6SJeremy L Thompson if (p < points_in_elem) { 443a2968d6SJeremy L Thompson const CeedInt ind = indices[p + elem * NUM_PTS]; 453a2968d6SJeremy L Thompson 463a2968d6SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 473a2968d6SJeremy L Thompson d_u[ind + comp * COMP_STRIDE] += r_u[comp]; 483a2968d6SJeremy L Thompson } 493a2968d6SJeremy L Thompson } 503a2968d6SJeremy L Thompson } 513a2968d6SJeremy L Thompson 523a2968d6SJeremy L Thompson //------------------------------------------------------------------------------ 539e201c85SYohann // 1D 549e201c85SYohann //------------------------------------------------------------------------------ 559e201c85SYohann 569e201c85SYohann //------------------------------------------------------------------------------ 579e201c85SYohann // L-vector -> E-vector, offsets provided 589e201c85SYohann //------------------------------------------------------------------------------ 596b92dc4bSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D> 60f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 61672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 626b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D) { 639e201c85SYohann const CeedInt node = data.t_id_x; 646b92dc4bSJeremy L Thompson const CeedInt ind = indices[node + elem * P_1D]; 65672b0f2aSSebastian Grimberg 66672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 679e201c85SYohann } 689e201c85SYohann } 699e201c85SYohann 709e201c85SYohann //------------------------------------------------------------------------------ 719e201c85SYohann // L-vector -> E-vector, strided 729e201c85SYohann //------------------------------------------------------------------------------ 736b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 74f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 756b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D) { 769e201c85SYohann const CeedInt node = data.t_id_x; 779e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 78672b0f2aSSebastian Grimberg 79672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 809e201c85SYohann } 819e201c85SYohann } 829e201c85SYohann 839e201c85SYohann //------------------------------------------------------------------------------ 849e201c85SYohann // E-vector -> L-vector, offsets provided 859e201c85SYohann //------------------------------------------------------------------------------ 866b92dc4bSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D> 87f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 88672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 896b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D) { 909e201c85SYohann const CeedInt node = data.t_id_x; 916b92dc4bSJeremy L Thompson const CeedInt ind = indices[node + elem * P_1D]; 92672b0f2aSSebastian Grimberg 93672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 949e201c85SYohann } 959e201c85SYohann } 969e201c85SYohann 979e201c85SYohann //------------------------------------------------------------------------------ 989e201c85SYohann // E-vector -> L-vector, strided 999e201c85SYohann //------------------------------------------------------------------------------ 1006b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 101f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 102672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 1036b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D) { 1049e201c85SYohann const CeedInt node = data.t_id_x; 1059e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 106672b0f2aSSebastian Grimberg 107672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 1089e201c85SYohann } 1099e201c85SYohann } 1109e201c85SYohann 1119e201c85SYohann //------------------------------------------------------------------------------ 1129e201c85SYohann // 2D 1139e201c85SYohann //------------------------------------------------------------------------------ 1149e201c85SYohann 1159e201c85SYohann //------------------------------------------------------------------------------ 1169e201c85SYohann // L-vector -> E-vector, offsets provided 1179e201c85SYohann //------------------------------------------------------------------------------ 1186b92dc4bSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D> 119f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 120672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 1216b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1226b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 1236b92dc4bSJeremy L Thompson const CeedInt ind = indices[node + elem * P_1D * P_1D]; 124672b0f2aSSebastian Grimberg 125672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 1269e201c85SYohann } 1279e201c85SYohann } 1289e201c85SYohann 1299e201c85SYohann //------------------------------------------------------------------------------ 1309e201c85SYohann // L-vector -> E-vector, strided 1319e201c85SYohann //------------------------------------------------------------------------------ 1326b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 133f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 1346b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1356b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 1369e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 137672b0f2aSSebastian Grimberg 138672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 1399e201c85SYohann } 1409e201c85SYohann } 1419e201c85SYohann 1429e201c85SYohann //------------------------------------------------------------------------------ 1439e201c85SYohann // E-vector -> L-vector, offsets provided 1449e201c85SYohann //------------------------------------------------------------------------------ 1456b92dc4bSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D> 146f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 147672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 1486b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1496b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 1506b92dc4bSJeremy L Thompson const CeedInt ind = indices[node + elem * P_1D * P_1D]; 151672b0f2aSSebastian Grimberg 152672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 1539e201c85SYohann } 1549e201c85SYohann } 1559e201c85SYohann 1569e201c85SYohann //------------------------------------------------------------------------------ 1579e201c85SYohann // E-vector -> L-vector, strided 1589e201c85SYohann //------------------------------------------------------------------------------ 1596b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 160f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 161672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 1626b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1636b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D; 1649e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 165672b0f2aSSebastian Grimberg 166672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 1679e201c85SYohann } 1689e201c85SYohann } 1699e201c85SYohann 1709e201c85SYohann //------------------------------------------------------------------------------ 1719e201c85SYohann // 3D 1729e201c85SYohann //------------------------------------------------------------------------------ 1739e201c85SYohann 1749e201c85SYohann //------------------------------------------------------------------------------ 1759e201c85SYohann // L-vector -> E-vector, offsets provided 1769e201c85SYohann //------------------------------------------------------------------------------ 1776b92dc4bSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D> 178f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 179672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 1806b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1816b92dc4bSJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 1826b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1836b92dc4bSJeremy L Thompson const CeedInt ind = indices[node + elem * P_1D * P_1D * P_1D]; 184672b0f2aSSebastian Grimberg 1856b92dc4bSJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1D] = d_u[ind + COMP_STRIDE * comp]; 1869e201c85SYohann } 1879e201c85SYohann } 188688b5473SJeremy L Thompson } 1899e201c85SYohann 1909e201c85SYohann //------------------------------------------------------------------------------ 1919e201c85SYohann // L-vector -> E-vector, strided 1929e201c85SYohann //------------------------------------------------------------------------------ 1936b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 194f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 1956b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 1966b92dc4bSJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 1976b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 1989e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 199672b0f2aSSebastian Grimberg 2006b92dc4bSJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1D] = d_u[ind + comp * STRIDES_COMP]; 2019e201c85SYohann } 2029e201c85SYohann } 203688b5473SJeremy L Thompson } 2049e201c85SYohann 2059e201c85SYohann //------------------------------------------------------------------------------ 2069e201c85SYohann // E-vector -> Q-vector, offests provided 2079e201c85SYohann //------------------------------------------------------------------------------ 2086b92dc4bSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int Q_1D> 209f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStandard3d(SharedData_Hip &data, const CeedInt nquads, const CeedInt elem, const CeedInt q, 210f815fac9SJeremy L Thompson const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, 211f815fac9SJeremy L Thompson CeedScalar *__restrict__ r_u) { 2126b92dc4bSJeremy L Thompson if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) { 2136b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * Q_1D + q * Q_1D * Q_1D; 2146b92dc4bSJeremy L Thompson const CeedInt ind = indices[node + elem * Q_1D * Q_1D * Q_1D]; 215672b0f2aSSebastian Grimberg 216672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 2179e201c85SYohann } 2189e201c85SYohann } 2199e201c85SYohann 2209e201c85SYohann //------------------------------------------------------------------------------ 2219e201c85SYohann // E-vector -> Q-vector, strided 2229e201c85SYohann //------------------------------------------------------------------------------ 2236b92dc4bSJeremy L Thompson template <int NUM_COMP, int Q_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 224f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u, 225672b0f2aSSebastian Grimberg CeedScalar *__restrict__ r_u) { 2266b92dc4bSJeremy L Thompson if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) { 2276b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * Q_1D + q * Q_1D * Q_1D; 2289e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 229672b0f2aSSebastian Grimberg 230672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 2319e201c85SYohann } 2329e201c85SYohann } 2339e201c85SYohann 2349e201c85SYohann //------------------------------------------------------------------------------ 2359e201c85SYohann // E-vector -> L-vector, offsets provided 2369e201c85SYohann //------------------------------------------------------------------------------ 2376b92dc4bSJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1D> 238f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 239672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 2406b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 2416b92dc4bSJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 2426b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 2436b92dc4bSJeremy L Thompson const CeedInt ind = indices[node + elem * P_1D * P_1D * P_1D]; 244672b0f2aSSebastian Grimberg 2456b92dc4bSJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1D]); 2469e201c85SYohann } 2479e201c85SYohann } 248688b5473SJeremy L Thompson } 2499e201c85SYohann 2509e201c85SYohann //------------------------------------------------------------------------------ 2519e201c85SYohann // E-vector -> L-vector, strided 2529e201c85SYohann //------------------------------------------------------------------------------ 2536b92dc4bSJeremy L Thompson template <int NUM_COMP, int P_1D, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 254f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 255672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 2566b92dc4bSJeremy L Thompson if (data.t_id_x < P_1D && data.t_id_y < P_1D) { 2576b92dc4bSJeremy L Thompson for (CeedInt z = 0; z < P_1D; z++) { 2586b92dc4bSJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1D + z * P_1D * P_1D; 2599e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 260672b0f2aSSebastian Grimberg 2616b92dc4bSJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1D]; 2629e201c85SYohann } 2639e201c85SYohann } 264688b5473SJeremy L Thompson } 2659e201c85SYohann 2669e201c85SYohann //------------------------------------------------------------------------------ 2679e201c85SYohann // 3D collocated derivatives computation 2689e201c85SYohann //------------------------------------------------------------------------------ 2696b92dc4bSJeremy L Thompson template <int NUM_COMP, int Q_1D, int T_1D> 270f815fac9SJeremy L Thompson inline __device__ void GradColloSlice3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 2712b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 2726b92dc4bSJeremy L Thompson if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) { 273672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 274*d6c19ee8SJeremy L Thompson __syncthreads(); 2756b92dc4bSJeremy L Thompson data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1D]; 2769e201c85SYohann __syncthreads(); 2779e201c85SYohann // X derivative 278672b0f2aSSebastian Grimberg r_V[comp + 0 * NUM_COMP] = 0.0; 2796b92dc4bSJeremy L Thompson for (CeedInt i = 0; i < Q_1D; i++) { 2806b92dc4bSJeremy L Thompson r_V[comp + 0 * NUM_COMP] += c_G[i + data.t_id_x * Q_1D] * data.slice[i + data.t_id_y * T_1D]; 281688b5473SJeremy L Thompson } 2829e201c85SYohann // Y derivative 283672b0f2aSSebastian Grimberg r_V[comp + 1 * NUM_COMP] = 0.0; 2846b92dc4bSJeremy L Thompson for (CeedInt i = 0; i < Q_1D; i++) { 2856b92dc4bSJeremy L Thompson r_V[comp + 1 * NUM_COMP] += c_G[i + data.t_id_y * Q_1D] * data.slice[data.t_id_x + i * T_1D]; 286688b5473SJeremy L Thompson } 2879e201c85SYohann // Z derivative 288672b0f2aSSebastian Grimberg r_V[comp + 2 * NUM_COMP] = 0.0; 2896b92dc4bSJeremy L Thompson for (CeedInt i = 0; i < Q_1D; i++) { 2906b92dc4bSJeremy L Thompson r_V[comp + 2 * NUM_COMP] += c_G[i + q * Q_1D] * r_U[i + comp * Q_1D]; 291688b5473SJeremy L Thompson } 2929e201c85SYohann } 2939e201c85SYohann } 2949e201c85SYohann } 2959e201c85SYohann 2969e201c85SYohann //------------------------------------------------------------------------------ 2979e201c85SYohann // 3D collocated derivatives transpose 2989e201c85SYohann //------------------------------------------------------------------------------ 2996b92dc4bSJeremy L Thompson template <int NUM_COMP, int Q_1D, int T_1D> 300f815fac9SJeremy L Thompson inline __device__ void GradColloSliceTranspose3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 3012b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 3026b92dc4bSJeremy L Thompson if (data.t_id_x < Q_1D && data.t_id_y < Q_1D) { 303672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 3049e201c85SYohann // X derivative 305*d6c19ee8SJeremy L Thompson __syncthreads(); 306672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP]; 3079e201c85SYohann __syncthreads(); 3086b92dc4bSJeremy L Thompson for (CeedInt i = 0; i < Q_1D; i++) { 3096b92dc4bSJeremy L Thompson r_V[q + comp * Q_1D] += c_G[data.t_id_x + i * Q_1D] * data.slice[i + data.t_id_y * T_1D]; 310688b5473SJeremy L Thompson } 3119e201c85SYohann // Y derivative 312*d6c19ee8SJeremy L Thompson __syncthreads(); 313672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP]; 3149e201c85SYohann __syncthreads(); 3156b92dc4bSJeremy L Thompson for (CeedInt i = 0; i < Q_1D; i++) { 3166b92dc4bSJeremy L Thompson r_V[q + comp * Q_1D] += c_G[data.t_id_y + i * Q_1D] * data.slice[data.t_id_x + i * T_1D]; 317688b5473SJeremy L Thompson } 3189e201c85SYohann // Z derivative 3196b92dc4bSJeremy L Thompson for (CeedInt i = 0; i < Q_1D; i++) { 3206b92dc4bSJeremy L Thompson r_V[i + comp * Q_1D] += c_G[i + q * Q_1D] * r_U[comp + 2 * NUM_COMP]; 321688b5473SJeremy L Thompson } 3229e201c85SYohann } 3239e201c85SYohann } 3249e201c85SYohann } 325