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 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> 16*f815fac9SJeremy 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 //------------------------------------------------------------------------------ 219e201c85SYohann // 1D 229e201c85SYohann //------------------------------------------------------------------------------ 239e201c85SYohann 249e201c85SYohann //------------------------------------------------------------------------------ 259e201c85SYohann // L-vector -> E-vector, offsets provided 269e201c85SYohann //------------------------------------------------------------------------------ 27672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 28*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 29672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 30672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 319e201c85SYohann const CeedInt node = data.t_id_x; 32672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d]; 33672b0f2aSSebastian Grimberg 34672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 359e201c85SYohann } 369e201c85SYohann } 379e201c85SYohann 389e201c85SYohann //------------------------------------------------------------------------------ 399e201c85SYohann // L-vector -> E-vector, strided 409e201c85SYohann //------------------------------------------------------------------------------ 41672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 42*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 43672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 449e201c85SYohann const CeedInt node = data.t_id_x; 459e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 46672b0f2aSSebastian Grimberg 47672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 489e201c85SYohann } 499e201c85SYohann } 509e201c85SYohann 519e201c85SYohann //------------------------------------------------------------------------------ 529e201c85SYohann // E-vector -> L-vector, offsets provided 539e201c85SYohann //------------------------------------------------------------------------------ 54672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 55*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard1d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 56672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 57672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 589e201c85SYohann const CeedInt node = data.t_id_x; 59672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d]; 60672b0f2aSSebastian Grimberg 61672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 629e201c85SYohann } 639e201c85SYohann } 649e201c85SYohann 659e201c85SYohann //------------------------------------------------------------------------------ 669e201c85SYohann // E-vector -> L-vector, strided 679e201c85SYohann //------------------------------------------------------------------------------ 68672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 69*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided1d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 70672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 71672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d) { 729e201c85SYohann const CeedInt node = data.t_id_x; 739e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 74672b0f2aSSebastian Grimberg 75672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 769e201c85SYohann } 779e201c85SYohann } 789e201c85SYohann 799e201c85SYohann //------------------------------------------------------------------------------ 809e201c85SYohann // 2D 819e201c85SYohann //------------------------------------------------------------------------------ 829e201c85SYohann 839e201c85SYohann //------------------------------------------------------------------------------ 849e201c85SYohann // L-vector -> E-vector, offsets provided 859e201c85SYohann //------------------------------------------------------------------------------ 86672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 87*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 88672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 89672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 90672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 91672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d]; 92672b0f2aSSebastian Grimberg 93672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 949e201c85SYohann } 959e201c85SYohann } 969e201c85SYohann 979e201c85SYohann //------------------------------------------------------------------------------ 989e201c85SYohann // L-vector -> E-vector, strided 999e201c85SYohann //------------------------------------------------------------------------------ 100672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 101*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 102672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 103672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 1049e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 105672b0f2aSSebastian Grimberg 106672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 1079e201c85SYohann } 1089e201c85SYohann } 1099e201c85SYohann 1109e201c85SYohann //------------------------------------------------------------------------------ 1119e201c85SYohann // E-vector -> L-vector, offsets provided 1129e201c85SYohann //------------------------------------------------------------------------------ 113672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 114*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard2d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 115672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 116672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 117672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 118672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d]; 119672b0f2aSSebastian Grimberg 120672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[comp]); 1219e201c85SYohann } 1229e201c85SYohann } 1239e201c85SYohann 1249e201c85SYohann //------------------------------------------------------------------------------ 1259e201c85SYohann // E-vector -> L-vector, strided 1269e201c85SYohann //------------------------------------------------------------------------------ 127672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 128*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided2d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 129672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 130672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) { 131672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d; 1329e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 133672b0f2aSSebastian Grimberg 134672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[comp]; 1359e201c85SYohann } 1369e201c85SYohann } 1379e201c85SYohann 1389e201c85SYohann //------------------------------------------------------------------------------ 1399e201c85SYohann // 3D 1409e201c85SYohann //------------------------------------------------------------------------------ 1419e201c85SYohann 1429e201c85SYohann //------------------------------------------------------------------------------ 1439e201c85SYohann // L-vector -> E-vector, offsets provided 1449e201c85SYohann //------------------------------------------------------------------------------ 145672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 146*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStandard3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 147672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 148672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 149672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 150672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 151672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 152672b0f2aSSebastian Grimberg 153672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + COMP_STRIDE * comp]; 1549e201c85SYohann } 1559e201c85SYohann } 1569e201c85SYohann 1579e201c85SYohann //------------------------------------------------------------------------------ 1589e201c85SYohann // L-vector -> E-vector, strided 1599e201c85SYohann //------------------------------------------------------------------------------ 160672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 161*f815fac9SJeremy L Thompson inline __device__ void ReadLVecStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *__restrict__ r_u) { 162672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 163672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 164672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 1659e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 166672b0f2aSSebastian Grimberg 167672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[z + comp * P_1d] = d_u[ind + comp * STRIDES_COMP]; 1689e201c85SYohann } 1699e201c85SYohann } 1709e201c85SYohann 1719e201c85SYohann //------------------------------------------------------------------------------ 1729e201c85SYohann // E-vector -> Q-vector, offests provided 1739e201c85SYohann //------------------------------------------------------------------------------ 174672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int Q_1d> 175*f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStandard3d(SharedData_Hip &data, const CeedInt nquads, const CeedInt elem, const CeedInt q, 176*f815fac9SJeremy L Thompson const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, 177*f815fac9SJeremy L Thompson CeedScalar *__restrict__ r_u) { 178672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 179672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 180672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * Q_1d * Q_1d * Q_1d]; 181672b0f2aSSebastian Grimberg 182672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + COMP_STRIDE * comp]; 1839e201c85SYohann } 1849e201c85SYohann } 1859e201c85SYohann 1869e201c85SYohann //------------------------------------------------------------------------------ 1879e201c85SYohann // E-vector -> Q-vector, strided 1889e201c85SYohann //------------------------------------------------------------------------------ 189672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 190*f815fac9SJeremy L Thompson inline __device__ void ReadEVecSliceStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u, 191672b0f2aSSebastian Grimberg CeedScalar *__restrict__ r_u) { 192672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 193672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * Q_1d + q * Q_1d * Q_1d; 1949e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 195672b0f2aSSebastian Grimberg 196672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) r_u[comp] = d_u[ind + comp * STRIDES_COMP]; 1979e201c85SYohann } 1989e201c85SYohann } 1999e201c85SYohann 2009e201c85SYohann //------------------------------------------------------------------------------ 2019e201c85SYohann // E-vector -> L-vector, offsets provided 2029e201c85SYohann //------------------------------------------------------------------------------ 203672b0f2aSSebastian Grimberg template <int NUM_COMP, int COMP_STRIDE, int P_1d> 204*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStandard3d(SharedData_Hip &data, const CeedInt num_nodes, const CeedInt elem, const CeedInt *__restrict__ indices, 205672b0f2aSSebastian Grimberg const CeedScalar *__restrict__ r_v, CeedScalar *__restrict__ d_v) { 206672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 207672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 208672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 209672b0f2aSSebastian Grimberg const CeedInt ind = indices[node + elem * P_1d * P_1d * P_1d]; 210672b0f2aSSebastian Grimberg 211672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) atomicAdd(&d_v[ind + COMP_STRIDE * comp], r_v[z + comp * P_1d]); 2129e201c85SYohann } 2139e201c85SYohann } 2149e201c85SYohann 2159e201c85SYohann //------------------------------------------------------------------------------ 2169e201c85SYohann // E-vector -> L-vector, strided 2179e201c85SYohann //------------------------------------------------------------------------------ 218672b0f2aSSebastian Grimberg template <int NUM_COMP, int P_1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM> 219*f815fac9SJeremy L Thompson inline __device__ void WriteLVecStrided3d(SharedData_Hip &data, const CeedInt elem, const CeedScalar *__restrict__ r_v, 220672b0f2aSSebastian Grimberg CeedScalar *__restrict__ d_v) { 221672b0f2aSSebastian Grimberg if (data.t_id_x < P_1d && data.t_id_y < P_1d) 222672b0f2aSSebastian Grimberg for (CeedInt z = 0; z < P_1d; z++) { 223672b0f2aSSebastian Grimberg const CeedInt node = data.t_id_x + data.t_id_y * P_1d + z * P_1d * P_1d; 2249e201c85SYohann const CeedInt ind = node * STRIDES_NODE + elem * STRIDES_ELEM; 225672b0f2aSSebastian Grimberg 226672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P_1d]; 2279e201c85SYohann } 2289e201c85SYohann } 2299e201c85SYohann 2309e201c85SYohann //------------------------------------------------------------------------------ 2319e201c85SYohann // 3D collocated derivatives computation 2329e201c85SYohann //------------------------------------------------------------------------------ 233672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d> 234*f815fac9SJeremy L Thompson inline __device__ void GradColloSlice3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 2352b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 236672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 237672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 238672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q_1d]; 2399e201c85SYohann __syncthreads(); 2409e201c85SYohann // X derivative 241672b0f2aSSebastian Grimberg r_V[comp + 0 * NUM_COMP] = 0.0; 242672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 243672b0f2aSSebastian 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) 2449e201c85SYohann // Y derivative 245672b0f2aSSebastian Grimberg r_V[comp + 1 * NUM_COMP] = 0.0; 246672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 247672b0f2aSSebastian 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) 2489e201c85SYohann // Z derivative 249672b0f2aSSebastian Grimberg r_V[comp + 2 * NUM_COMP] = 0.0; 250672b0f2aSSebastian 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) 2519e201c85SYohann __syncthreads(); 2529e201c85SYohann } 2539e201c85SYohann } 2549e201c85SYohann } 2559e201c85SYohann 2569e201c85SYohann //------------------------------------------------------------------------------ 2579e201c85SYohann // 3D collocated derivatives transpose 2589e201c85SYohann //------------------------------------------------------------------------------ 259672b0f2aSSebastian Grimberg template <int NUM_COMP, int Q_1d> 260*f815fac9SJeremy L Thompson inline __device__ void GradColloSliceTranspose3d(SharedData_Hip &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, 2612b730f8bSJeremy L Thompson CeedScalar *__restrict__ r_V) { 262672b0f2aSSebastian Grimberg if (data.t_id_x < Q_1d && data.t_id_y < Q_1d) { 263672b0f2aSSebastian Grimberg for (CeedInt comp = 0; comp < NUM_COMP; comp++) { 2649e201c85SYohann // X derivative 265672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NUM_COMP]; 2669e201c85SYohann __syncthreads(); 267672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 268672b0f2aSSebastian 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) 2699e201c85SYohann __syncthreads(); 2709e201c85SYohann // Y derivative 271672b0f2aSSebastian Grimberg data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NUM_COMP]; 2729e201c85SYohann __syncthreads(); 273672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 274672b0f2aSSebastian 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) 2759e201c85SYohann __syncthreads(); 2769e201c85SYohann // Z derivative 277672b0f2aSSebastian Grimberg for (CeedInt i = 0; i < Q_1d; i++) 278672b0f2aSSebastian 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) 2799e201c85SYohann } 2809e201c85SYohann } 2819e201c85SYohann } 282