15aed82e4SJeremy 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 CUDA backend macro and type definitions for JiT source 109e201c85SYohann 11*c0b5abf0SJeremy L Thompson #include <ceed/types.h> 129e201c85SYohann 139e201c85SYohann //------------------------------------------------------------------------------ 149e201c85SYohann // Load matrices for basis actions 159e201c85SYohann //------------------------------------------------------------------------------ 169e201c85SYohann template <int P, int Q> 179e201c85SYohann inline __device__ void loadMatrix(SharedData_Cuda &data, const CeedScalar *__restrict__ d_B, CeedScalar *B) { 182b730f8bSJeremy L Thompson for (CeedInt i = data.t_id; i < P * Q; i += blockDim.x * blockDim.y * blockDim.z) B[i] = d_B[i]; 199e201c85SYohann } 209e201c85SYohann 219e201c85SYohann //------------------------------------------------------------------------------ 229e201c85SYohann // 1D 239e201c85SYohann //------------------------------------------------------------------------------ 249e201c85SYohann 259e201c85SYohann //------------------------------------------------------------------------------ 269e201c85SYohann // L-vector -> E-vector, offsets provided 279e201c85SYohann //------------------------------------------------------------------------------ 28ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d> 29ca735530SJeremy L Thompson inline __device__ void readDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 30672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 31ca735530SJeremy L Thompson if (data.t_id_x < P_1d) { 329e201c85SYohann const CeedInt node = data.t_id_x; 33ca735530SJeremy L Thompson const CeedInt ind = indices[node + elem * P_1d]; 34ca735530SJeremy L Thompson 35ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 369e201c85SYohann } 379e201c85SYohann } 389e201c85SYohann 399e201c85SYohann //------------------------------------------------------------------------------ 409e201c85SYohann // L-vector -> E-vector, strided 419e201c85SYohann //------------------------------------------------------------------------------ 42ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 43672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, 44672b0f2aSSebastian Grimberg CeedScalar *__restrict__ r_u) { 45ca735530SJeremy L Thompson if (data.t_id_x < P_1d) { 469e201c85SYohann const CeedInt node = data.t_id_x; 479e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 48ca735530SJeremy L Thompson 49ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 509e201c85SYohann } 519e201c85SYohann } 529e201c85SYohann 539e201c85SYohann //------------------------------------------------------------------------------ 549e201c85SYohann // E-vector -> L-vector, offsets provided 559e201c85SYohann //------------------------------------------------------------------------------ 56ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d> 57ca735530SJeremy L Thompson inline __device__ void writeDofsOffset1d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 58672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 59ca735530SJeremy L Thompson if (data.t_id_x < P_1d) { 609e201c85SYohann const CeedInt node = data.t_id_x; 61ca735530SJeremy L Thompson const CeedInt ind = indices[node + elem * P_1d]; 62ca735530SJeremy L Thompson 63ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 649e201c85SYohann } 659e201c85SYohann } 669e201c85SYohann 679e201c85SYohann //------------------------------------------------------------------------------ 689e201c85SYohann // E-vector -> L-vector, strided 699e201c85SYohann //------------------------------------------------------------------------------ 70ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 71672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 72672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 73ca735530SJeremy L Thompson if (data.t_id_x < P_1d) { 749e201c85SYohann const CeedInt node = data.t_id_x; 759e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 76ca735530SJeremy L Thompson 77ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 789e201c85SYohann } 799e201c85SYohann } 809e201c85SYohann 819e201c85SYohann //------------------------------------------------------------------------------ 829e201c85SYohann // 2D 839e201c85SYohann //------------------------------------------------------------------------------ 849e201c85SYohann 859e201c85SYohann //------------------------------------------------------------------------------ 869e201c85SYohann // L-vector -> E-vector, offsets provided 879e201c85SYohann //------------------------------------------------------------------------------ 88ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d> 89ca735530SJeremy L Thompson inline __device__ void readDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 90672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 91ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 92ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 93ca735530SJeremy L Thompson const CeedInt ind = indices[node + elem * P_1d * P_1d]; 94ca735530SJeremy L Thompson 95ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 969e201c85SYohann } 979e201c85SYohann } 989e201c85SYohann 999e201c85SYohann //------------------------------------------------------------------------------ 1009e201c85SYohann // L-vector -> E-vector, strided 1019e201c85SYohann //------------------------------------------------------------------------------ 102ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 103672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, 104672b0f2aSSebastian Grimberg CeedScalar *__restrict__ r_u) { 105ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 106ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 1079e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 108ca735530SJeremy L Thompson 109ca735530SJeremy L Thompson 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 //------------------------------------------------------------------------------ 116ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d> 117ca735530SJeremy L Thompson inline __device__ void writeDofsOffset2d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 118672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 119ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 120ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 121ca735530SJeremy L Thompson const CeedInt ind = indices[node + elem * P_1d * P_1d]; 122ca735530SJeremy L Thompson 123ca735530SJeremy L Thompson 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 //------------------------------------------------------------------------------ 130ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 131672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 132672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 133ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 134ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 1359e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 136ca735530SJeremy L Thompson 137ca735530SJeremy L Thompson 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 //------------------------------------------------------------------------------ 1489e201c85SYohann // TODO: remove "Dofs" and "Quads" in the following function names? 1499e201c85SYohann // - readDofsOffset3d -> readOffset3d ? 1509e201c85SYohann // - readDofsStrided3d -> readStrided3d ? 1519e201c85SYohann // - readSliceQuadsOffset3d -> readSliceOffset3d ? 1529e201c85SYohann // - readSliceQuadsStrided3d -> readSliceStrided3d ? 1539e201c85SYohann // - writeDofsOffset3d -> writeOffset3d ? 1549e201c85SYohann // - writeDofsStrided3d -> writeStrided3d ? 155ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d> 156ca735530SJeremy L Thompson inline __device__ void readDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 157672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 158ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) 159ca735530SJeremy L Thompson for (CeedInt z = 0; z < P_1d; z++) { 160ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 161ca735530SJeremy L Thompson const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 162ca735530SJeremy L Thompson 163ca735530SJeremy L Thompson 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 //------------------------------------------------------------------------------ 170ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 171672b0f2aSSebastian Grimberg inline __device__ void readDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, 172672b0f2aSSebastian Grimberg CeedScalar *__restrict__ r_u) { 173ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) 174ca735530SJeremy L Thompson for (CeedInt z = 0; z < P_1d; z++) { 175ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 1769e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 177ca735530SJeremy L Thompson 178ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP]; 1799e201c85SYohann } 1809e201c85SYohann } 1819e201c85SYohann 1829e201c85SYohann //------------------------------------------------------------------------------ 1839e201c85SYohann // E-vector -> Q-vector, offests provided 1849e201c85SYohann //------------------------------------------------------------------------------ 185ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int Q_1d> 1862b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsOffset3d(SharedData_Cuda &data, const CeedInt nquads, const CeedInt elem, const CeedInt q, 187672b0f2aSSebastian Grimberg const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 188ca735530SJeremy L Thompson if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 189ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 190ca735530SJeremy L Thompson const CeedInt ind = indices[node + elem * Q_1d * Q_1d * Q_1d]; 191ca735530SJeremy L Thompson 192ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 1939e201c85SYohann } 1949e201c85SYohann } 1959e201c85SYohann 1969e201c85SYohann //------------------------------------------------------------------------------ 1979e201c85SYohann // E-vector -> Q-vector, strided 1989e201c85SYohann //------------------------------------------------------------------------------ 199ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 2002b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u, 201672b0f2aSSebastian Grimberg CeedScalar *__restrict__ r_u) { 202ca735530SJeremy L Thompson if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 203ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 2049e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 205ca735530SJeremy L Thompson 206ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 2079e201c85SYohann } 2089e201c85SYohann } 2099e201c85SYohann 2109e201c85SYohann //------------------------------------------------------------------------------ 2119e201c85SYohann // E-vector -> L-vector, offsets provided 2129e201c85SYohann //------------------------------------------------------------------------------ 213ca735530SJeremy L Thompson template <int NUM_COMP, int COMP_STRIDE, int P_1d> 214ca735530SJeremy L Thompson inline __device__ void writeDofsOffset3d(SharedData_Cuda &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 215672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 216ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) 217ca735530SJeremy L Thompson for (CeedInt z = 0; z < P_1d; z++) { 218ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 219ca735530SJeremy L Thompson const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 220ca735530SJeremy L Thompson 221ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]); 2229e201c85SYohann } 2239e201c85SYohann } 2249e201c85SYohann 2259e201c85SYohann //------------------------------------------------------------------------------ 2269e201c85SYohann // E-vector -> L-vector, strided 2279e201c85SYohann //------------------------------------------------------------------------------ 228ca735530SJeremy L Thompson template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 229672b0f2aSSebastian Grimberg inline __device__ void writeDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 230672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 231ca735530SJeremy L Thompson if (data.t_id_x < P_1d && data.t_id_y < P_1d) 232ca735530SJeremy L Thompson for (CeedInt z = 0; z < P_1d; z++) { 233ca735530SJeremy L Thompson const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 2349e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 235ca735530SJeremy L Thompson 236ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d]; 2379e201c85SYohann } 2389e201c85SYohann } 2399e201c85SYohann 2409e201c85SYohann //------------------------------------------------------------------------------ 2419e201c85SYohann // 3D collocated derivatives computation 2429e201c85SYohann //------------------------------------------------------------------------------ 243ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d> 2442b730f8bSJeremy L Thompson inline __device__ void gradCollo3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 2452b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 246ca735530SJeremy L Thompson if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 247ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 248ca735530SJeremy L Thompson data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d]; 2499e201c85SYohann __syncthreads(); 2509e201c85SYohann // X derivative 251ca735530SJeremy L Thompson r_V[comp + 0 * NUM_COMP] = 0.0; 252ca735530SJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) 253ca735530SJeremy 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]; // Contract x direction (X derivative) 2549e201c85SYohann // Y derivative 255ca735530SJeremy L Thompson r_V[comp + 1 * NUM_COMP] = 0.0; 256ca735530SJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) 257ca735530SJeremy 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]; // Contract y direction (Y derivative) 2589e201c85SYohann // Z derivative 259ca735530SJeremy L Thompson r_V[comp + 2 * NUM_COMP] = 0.0; 260ca735530SJeremy L Thompson 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) 2619e201c85SYohann __syncthreads(); 2629e201c85SYohann } 2639e201c85SYohann } 2649e201c85SYohann } 2659e201c85SYohann 2669e201c85SYohann //------------------------------------------------------------------------------ 2679e201c85SYohann // 3D collocated derivatives transpose 2689e201c85SYohann //------------------------------------------------------------------------------ 269ca735530SJeremy L Thompson template <int NUM_COMP, int Q_1d> 2702b730f8bSJeremy L Thompson inline __device__ void gradColloTranspose3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 2712b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 272ca735530SJeremy L Thompson if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 273ca735530SJeremy L Thompson for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 2749e201c85SYohann // X derivative 275ca735530SJeremy L Thompson data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP]; 2769e201c85SYohann __syncthreads(); 277ca735530SJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) 278ca735530SJeremy 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]; // Contract x direction (X derivative) 2799e201c85SYohann __syncthreads(); 2809e201c85SYohann // Y derivative 281ca735530SJeremy L Thompson data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP]; 2829e201c85SYohann __syncthreads(); 283ca735530SJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) 284ca735530SJeremy 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]; // Contract y direction (Y derivative) 2859e201c85SYohann __syncthreads(); 2869e201c85SYohann // Z derivative 287ca735530SJeremy L Thompson for (CeedInt i = 0; i < Q_1d; i++) 288ca735530SJeremy L Thompson r_V[i + comp * Q_1d] += c_G[i + q * Q_1d] * r_U[comp + 2 * NUM_COMP]; // PARTIAL contract z direction (Z derivative) 2899e201c85SYohann } 2909e201c85SYohann } 2919e201c85SYohann } 292