xref: /libCEED/rust/libceed-sys/c-src/include/ceed/jit-source/cuda/cuda-gen-templates.h (revision 94b7b29b41ad8a17add4c577886859ef16f89dec)
19e201c85SYohann // Copyright (c) 2017-2022, 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
10*94b7b29bSJeremy L Thompson #ifndef CEED_CUDA_GEN_TEMPLATES_H
11*94b7b29bSJeremy L Thompson #define CEED_CUDA_GEN_TEMPLATES_H
129e201c85SYohann 
139e201c85SYohann #include <ceed/types.h>
149e201c85SYohann 
159e201c85SYohann //------------------------------------------------------------------------------
169e201c85SYohann // Load matrices for basis actions
179e201c85SYohann //------------------------------------------------------------------------------
189e201c85SYohann template <int P, int Q>
199e201c85SYohann inline __device__ void loadMatrix(SharedData_Cuda &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 //------------------------------------------------------------------------------
309e201c85SYohann template <int NCOMP, int COMPSTRIDE, int P1d>
312b730f8bSJeremy L Thompson inline __device__ void readDofsOffset1d(SharedData_Cuda &data, const CeedInt nnodes, const CeedInt elem, const CeedInt *__restrict__ indices,
322b730f8bSJeremy L Thompson                                         const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
339e201c85SYohann   if (data.t_id_x < P1d) {
349e201c85SYohann     const CeedInt node = data.t_id_x;
359e201c85SYohann     const CeedInt ind  = indices[node + elem * P1d];
362b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[comp] = d_u[ind + COMPSTRIDE * comp];
379e201c85SYohann   }
389e201c85SYohann }
399e201c85SYohann 
409e201c85SYohann //------------------------------------------------------------------------------
419e201c85SYohann // L-vector -> E-vector, strided
429e201c85SYohann //------------------------------------------------------------------------------
439e201c85SYohann template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
449e201c85SYohann inline __device__ void readDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
459e201c85SYohann   if (data.t_id_x < P1d) {
469e201c85SYohann     const CeedInt node = data.t_id_x;
479e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
482b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
499e201c85SYohann   }
509e201c85SYohann }
519e201c85SYohann 
529e201c85SYohann //------------------------------------------------------------------------------
539e201c85SYohann // E-vector -> L-vector, offsets provided
549e201c85SYohann //------------------------------------------------------------------------------
559e201c85SYohann template <int NCOMP, int COMPSTRIDE, int P1d>
562b730f8bSJeremy L Thompson inline __device__ void writeDofsOffset1d(SharedData_Cuda &data, const CeedInt nnodes, const CeedInt elem, const CeedInt *__restrict__ indices,
572b730f8bSJeremy L Thompson                                          const CeedScalar *r_v, CeedScalar *d_v) {
589e201c85SYohann   if (data.t_id_x < P1d) {
599e201c85SYohann     const CeedInt node = data.t_id_x;
609e201c85SYohann     const CeedInt ind  = indices[node + elem * P1d];
612b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) atomicAdd(&d_v[ind + COMPSTRIDE * comp], r_v[comp]);
629e201c85SYohann   }
639e201c85SYohann }
649e201c85SYohann 
659e201c85SYohann //------------------------------------------------------------------------------
669e201c85SYohann // E-vector -> L-vector, strided
679e201c85SYohann //------------------------------------------------------------------------------
689e201c85SYohann template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
699e201c85SYohann inline __device__ void writeDofsStrided1d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *r_v, CeedScalar *d_v) {
709e201c85SYohann   if (data.t_id_x < P1d) {
719e201c85SYohann     const CeedInt node = data.t_id_x;
729e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
732b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
749e201c85SYohann   }
759e201c85SYohann }
769e201c85SYohann 
779e201c85SYohann //------------------------------------------------------------------------------
789e201c85SYohann // 2D
799e201c85SYohann //------------------------------------------------------------------------------
809e201c85SYohann 
819e201c85SYohann //------------------------------------------------------------------------------
829e201c85SYohann // L-vector -> E-vector, offsets provided
839e201c85SYohann //------------------------------------------------------------------------------
849e201c85SYohann template <int NCOMP, int COMPSTRIDE, int P1d>
852b730f8bSJeremy L Thompson inline __device__ void readDofsOffset2d(SharedData_Cuda &data, const CeedInt nnodes, const CeedInt elem, const CeedInt *__restrict__ indices,
862b730f8bSJeremy L Thompson                                         const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
879e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d) {
889e201c85SYohann     const CeedInt node = data.t_id_x + data.t_id_y * P1d;
899e201c85SYohann     const CeedInt ind  = indices[node + elem * P1d * P1d];
902b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[comp] = d_u[ind + COMPSTRIDE * comp];
919e201c85SYohann   }
929e201c85SYohann }
939e201c85SYohann 
949e201c85SYohann //------------------------------------------------------------------------------
959e201c85SYohann // L-vector -> E-vector, strided
969e201c85SYohann //------------------------------------------------------------------------------
979e201c85SYohann template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
989e201c85SYohann inline __device__ void readDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
999e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d) {
1009e201c85SYohann     const CeedInt node = data.t_id_x + data.t_id_y * P1d;
1019e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
1022b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
1039e201c85SYohann   }
1049e201c85SYohann }
1059e201c85SYohann 
1069e201c85SYohann //------------------------------------------------------------------------------
1079e201c85SYohann // E-vector -> L-vector, offsets provided
1089e201c85SYohann //------------------------------------------------------------------------------
1099e201c85SYohann template <int NCOMP, int COMPSTRIDE, int P1d>
1102b730f8bSJeremy L Thompson inline __device__ void writeDofsOffset2d(SharedData_Cuda &data, const CeedInt nnodes, const CeedInt elem, const CeedInt *__restrict__ indices,
1112b730f8bSJeremy L Thompson                                          const CeedScalar *r_v, CeedScalar *d_v) {
1129e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d) {
1139e201c85SYohann     const CeedInt node = data.t_id_x + data.t_id_y * P1d;
1149e201c85SYohann     const CeedInt ind  = indices[node + elem * P1d * P1d];
1152b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) atomicAdd(&d_v[ind + COMPSTRIDE * comp], r_v[comp]);
1169e201c85SYohann   }
1179e201c85SYohann }
1189e201c85SYohann 
1199e201c85SYohann //------------------------------------------------------------------------------
1209e201c85SYohann // E-vector -> L-vector, strided
1219e201c85SYohann //------------------------------------------------------------------------------
1229e201c85SYohann template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
1239e201c85SYohann inline __device__ void writeDofsStrided2d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *r_v, CeedScalar *d_v) {
1249e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d) {
1259e201c85SYohann     const CeedInt node = data.t_id_x + data.t_id_y * P1d;
1269e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
1272b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) d_v[ind + comp * STRIDES_COMP] += r_v[comp];
1289e201c85SYohann   }
1299e201c85SYohann }
1309e201c85SYohann 
1319e201c85SYohann //------------------------------------------------------------------------------
1329e201c85SYohann // 3D
1339e201c85SYohann //------------------------------------------------------------------------------
1349e201c85SYohann 
1359e201c85SYohann //------------------------------------------------------------------------------
1369e201c85SYohann // L-vector -> E-vector, offsets provided
1379e201c85SYohann //------------------------------------------------------------------------------
1389e201c85SYohann // TODO: remove "Dofs" and "Quads" in the following function names?
1399e201c85SYohann //   - readDofsOffset3d -> readOffset3d ?
1409e201c85SYohann //   - readDofsStrided3d -> readStrided3d ?
1419e201c85SYohann //   - readSliceQuadsOffset3d -> readSliceOffset3d ?
1429e201c85SYohann //   - readSliceQuadsStrided3d -> readSliceStrided3d ?
1439e201c85SYohann //   - writeDofsOffset3d -> writeOffset3d ?
1449e201c85SYohann //   - writeDofsStrided3d -> writeStrided3d ?
1459e201c85SYohann template <int NCOMP, int COMPSTRIDE, int P1d>
1462b730f8bSJeremy L Thompson inline __device__ void readDofsOffset3d(SharedData_Cuda &data, const CeedInt nnodes, const CeedInt elem, const CeedInt *__restrict__ indices,
1472b730f8bSJeremy L Thompson                                         const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
1489e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d)
1499e201c85SYohann     for (CeedInt z = 0; z < P1d; ++z) {
1509e201c85SYohann       const CeedInt node = data.t_id_x + data.t_id_y * P1d + z * P1d * P1d;
1519e201c85SYohann       const CeedInt ind  = indices[node + elem * P1d * P1d * P1d];
1522b730f8bSJeremy L Thompson       for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[z + comp * P1d] = d_u[ind + COMPSTRIDE * comp];
1539e201c85SYohann     }
1549e201c85SYohann }
1559e201c85SYohann 
1569e201c85SYohann //------------------------------------------------------------------------------
1579e201c85SYohann // L-vector -> E-vector, strided
1589e201c85SYohann //------------------------------------------------------------------------------
1599e201c85SYohann template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
1609e201c85SYohann inline __device__ void readDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
1619e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d)
1629e201c85SYohann     for (CeedInt z = 0; z < P1d; ++z) {
1639e201c85SYohann       const CeedInt node = data.t_id_x + data.t_id_y * P1d + z * P1d * P1d;
1649e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
1652b730f8bSJeremy L Thompson       for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[z + comp * P1d] = d_u[ind + comp * STRIDES_COMP];
1669e201c85SYohann     }
1679e201c85SYohann }
1689e201c85SYohann 
1699e201c85SYohann //------------------------------------------------------------------------------
1709e201c85SYohann // E-vector -> Q-vector, offests provided
1719e201c85SYohann //------------------------------------------------------------------------------
1729e201c85SYohann template <int NCOMP, int COMPSTRIDE, int Q1d>
1732b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsOffset3d(SharedData_Cuda &data, const CeedInt nquads, const CeedInt elem, const CeedInt q,
1742b730f8bSJeremy L Thompson                                               const CeedInt *__restrict__ indices, const CeedScalar *__restrict__ d_u, CeedScalar *r_u) {
1759e201c85SYohann   if (data.t_id_x < Q1d && data.t_id_y < Q1d) {
1769e201c85SYohann     const CeedInt node = data.t_id_x + data.t_id_y * Q1d + q * Q1d * Q1d;
1772b730f8bSJeremy L Thompson     const CeedInt ind  = indices[node + elem * Q1d * Q1d * Q1d];
1782b730f8bSJeremy L Thompson     ;
1792b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[comp] = d_u[ind + COMPSTRIDE * comp];
1809e201c85SYohann   }
1819e201c85SYohann }
1829e201c85SYohann 
1839e201c85SYohann //------------------------------------------------------------------------------
1849e201c85SYohann // E-vector -> Q-vector, strided
1859e201c85SYohann //------------------------------------------------------------------------------
1869e201c85SYohann template <int NCOMP, int Q1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
1872b730f8bSJeremy L Thompson inline __device__ void readSliceQuadsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedInt q, const CeedScalar *__restrict__ d_u,
1882b730f8bSJeremy L Thompson                                                CeedScalar *r_u) {
1899e201c85SYohann   if (data.t_id_x < Q1d && data.t_id_y < Q1d) {
1909e201c85SYohann     const CeedInt node = data.t_id_x + data.t_id_y * Q1d + q * Q1d * Q1d;
1919e201c85SYohann     const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
1922b730f8bSJeremy L Thompson     for (CeedInt comp = 0; comp < NCOMP; ++comp) r_u[comp] = d_u[ind + comp * STRIDES_COMP];
1939e201c85SYohann   }
1949e201c85SYohann }
1959e201c85SYohann 
1969e201c85SYohann //------------------------------------------------------------------------------
1979e201c85SYohann // E-vector -> L-vector, offsets provided
1989e201c85SYohann //------------------------------------------------------------------------------
1999e201c85SYohann template <int NCOMP, int COMPSTRIDE, int P1d>
2002b730f8bSJeremy L Thompson inline __device__ void writeDofsOffset3d(SharedData_Cuda &data, const CeedInt nnodes, const CeedInt elem, const CeedInt *__restrict__ indices,
2012b730f8bSJeremy L Thompson                                          const CeedScalar *r_v, CeedScalar *d_v) {
2029e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d)
2039e201c85SYohann     for (CeedInt z = 0; z < P1d; ++z) {
2049e201c85SYohann       const CeedInt node = data.t_id_x + data.t_id_y * P1d + z * P1d * P1d;
2059e201c85SYohann       const CeedInt ind  = indices[node + elem * P1d * P1d * P1d];
2062b730f8bSJeremy L Thompson       for (CeedInt comp = 0; comp < NCOMP; ++comp) atomicAdd(&d_v[ind + COMPSTRIDE * comp], r_v[z + comp * P1d]);
2079e201c85SYohann     }
2089e201c85SYohann }
2099e201c85SYohann 
2109e201c85SYohann //------------------------------------------------------------------------------
2119e201c85SYohann // E-vector -> L-vector, strided
2129e201c85SYohann //------------------------------------------------------------------------------
2139e201c85SYohann template <int NCOMP, int P1d, int STRIDES_NODE, int STRIDES_COMP, int STRIDES_ELEM>
2149e201c85SYohann inline __device__ void writeDofsStrided3d(SharedData_Cuda &data, const CeedInt elem, const CeedScalar *r_v, CeedScalar *d_v) {
2159e201c85SYohann   if (data.t_id_x < P1d && data.t_id_y < P1d)
2169e201c85SYohann     for (CeedInt z = 0; z < P1d; ++z) {
2179e201c85SYohann       const CeedInt node = data.t_id_x + data.t_id_y * P1d + z * P1d * P1d;
2189e201c85SYohann       const CeedInt ind  = node * STRIDES_NODE + elem * STRIDES_ELEM;
2192b730f8bSJeremy L Thompson       for (CeedInt comp = 0; comp < NCOMP; ++comp) d_v[ind + comp * STRIDES_COMP] += r_v[z + comp * P1d];
2209e201c85SYohann     }
2219e201c85SYohann }
2229e201c85SYohann 
2239e201c85SYohann //------------------------------------------------------------------------------
2249e201c85SYohann // 3D collocated derivatives computation
2259e201c85SYohann //------------------------------------------------------------------------------
2269e201c85SYohann template <int NCOMP, int Q1d>
2272b730f8bSJeremy L Thompson inline __device__ void gradCollo3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2282b730f8bSJeremy L Thompson                                    CeedScalar *__restrict__ r_V) {
2299e201c85SYohann   if (data.t_id_x < Q1d && data.t_id_y < Q1d) {
2309e201c85SYohann     for (CeedInt comp = 0; comp < NCOMP; ++comp) {
2319e201c85SYohann       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[q + comp * Q1d];
2329e201c85SYohann       __syncthreads();
2339e201c85SYohann       // X derivative
2349e201c85SYohann       r_V[comp + 0 * NCOMP] = 0.0;
2359e201c85SYohann       for (CeedInt i = 0; i < Q1d; ++i)
2369e201c85SYohann         r_V[comp + 0 * NCOMP] += c_G[i + data.t_id_x * Q1d] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction (X derivative)
2379e201c85SYohann       // Y derivative
2389e201c85SYohann       r_V[comp + 1 * NCOMP] = 0.0;
2399e201c85SYohann       for (CeedInt i = 0; i < Q1d; ++i)
2409e201c85SYohann         r_V[comp + 1 * NCOMP] += c_G[i + data.t_id_y * Q1d] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction (Y derivative)
2419e201c85SYohann       // Z derivative
2429e201c85SYohann       r_V[comp + 2 * NCOMP] = 0.0;
2432b730f8bSJeremy L Thompson       for (CeedInt i = 0; i < Q1d; ++i) r_V[comp + 2 * NCOMP] += c_G[i + q * Q1d] * r_U[i + comp * Q1d];  // Contract z direction (Z derivative)
2449e201c85SYohann       __syncthreads();
2459e201c85SYohann     }
2469e201c85SYohann   }
2479e201c85SYohann }
2489e201c85SYohann 
2499e201c85SYohann //------------------------------------------------------------------------------
2509e201c85SYohann // 3D collocated derivatives transpose
2519e201c85SYohann //------------------------------------------------------------------------------
2529e201c85SYohann template <int NCOMP, int Q1d>
2532b730f8bSJeremy L Thompson inline __device__ void gradColloTranspose3d(SharedData_Cuda &data, const CeedInt q, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
2542b730f8bSJeremy L Thompson                                             CeedScalar *__restrict__ r_V) {
2559e201c85SYohann   if (data.t_id_x < Q1d && data.t_id_y < Q1d) {
2569e201c85SYohann     for (CeedInt comp = 0; comp < NCOMP; ++comp) {
2579e201c85SYohann       // X derivative
2589e201c85SYohann       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 0 * NCOMP];
2599e201c85SYohann       __syncthreads();
2609e201c85SYohann       for (CeedInt i = 0; i < Q1d; ++i)
2619e201c85SYohann         r_V[q + comp * Q1d] += c_G[data.t_id_x + i * Q1d] * data.slice[i + data.t_id_y * T_1D];  // Contract x direction (X derivative)
2629e201c85SYohann       __syncthreads();
2639e201c85SYohann       // Y derivative
2649e201c85SYohann       data.slice[data.t_id_x + data.t_id_y * T_1D] = r_U[comp + 1 * NCOMP];
2659e201c85SYohann       __syncthreads();
2669e201c85SYohann       for (CeedInt i = 0; i < Q1d; ++i)
2679e201c85SYohann         r_V[q + comp * Q1d] += c_G[data.t_id_y + i * Q1d] * data.slice[data.t_id_x + i * T_1D];  // Contract y direction (Y derivative)
2689e201c85SYohann       __syncthreads();
2699e201c85SYohann       // Z derivative
2709e201c85SYohann       for (CeedInt i = 0; i < Q1d; ++i)
2719e201c85SYohann         r_V[i + comp * Q1d] += c_G[i + q * Q1d] * r_U[comp + 2 * NCOMP];  // PARTIAL contract z direction (Z derivative)
2729e201c85SYohann     }
2739e201c85SYohann   }
2749e201c85SYohann }
2759e201c85SYohann 
276*94b7b29bSJeremy L Thompson #endif  // CEED_CUDA_GEN_TEMPLATES_H
277