xref: /libCEED/rust/libceed-sys/c-src/backends/hip-ref/ceed-hip-ref-restriction.c (revision 6574a04ff2135c3834f1b6ef9a4ec7566c4782db)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
849aac155SJeremy L Thompson #include <ceed.h>
90d0321e0SJeremy L Thompson #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
110d0321e0SJeremy L Thompson #include <hip/hip_runtime.h>
120d0321e0SJeremy L Thompson #include <stdbool.h>
130d0321e0SJeremy L Thompson #include <stddef.h>
1444d7a66cSJeremy L Thompson #include <string.h>
152b730f8bSJeremy L Thompson 
1649aac155SJeremy L Thompson #include "../hip/ceed-hip-common.h"
170d0321e0SJeremy L Thompson #include "../hip/ceed-hip-compile.h"
182b730f8bSJeremy L Thompson #include "ceed-hip-ref.h"
190d0321e0SJeremy L Thompson 
200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
210d0321e0SJeremy L Thompson // Apply restriction
220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
232b730f8bSJeremy L Thompson static int CeedElemRestrictionApply_Hip(CeedElemRestriction r, CeedTransposeMode t_mode, CeedVector u, CeedVector v, CeedRequest *request) {
240d0321e0SJeremy L Thompson   CeedElemRestriction_Hip *impl;
252b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(r, &impl));
260d0321e0SJeremy L Thompson   Ceed ceed;
272b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
280d0321e0SJeremy L Thompson   Ceed_Hip *data;
292b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &data));
30437930d1SJeremy L Thompson   const CeedInt block_size = 64;
31437930d1SJeremy L Thompson   const CeedInt num_nodes  = impl->num_nodes;
32437930d1SJeremy L Thompson   CeedInt       num_elem, elem_size;
33437930d1SJeremy L Thompson   CeedElemRestrictionGetNumElements(r, &num_elem);
342b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size));
350d0321e0SJeremy L Thompson   hipFunction_t kernel;
360d0321e0SJeremy L Thompson 
370d0321e0SJeremy L Thompson   // Get vectors
380d0321e0SJeremy L Thompson   const CeedScalar *d_u;
390d0321e0SJeremy L Thompson   CeedScalar       *d_v;
402b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
41437930d1SJeremy L Thompson   if (t_mode == CEED_TRANSPOSE) {
420d0321e0SJeremy L Thompson     // Sum into for transpose mode, e-vec to l-vec
432b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
440d0321e0SJeremy L Thompson   } else {
450d0321e0SJeremy L Thompson     // Overwrite for notranspose mode, l-vec to e-vec
462b730f8bSJeremy L Thompson     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
470d0321e0SJeremy L Thompson   }
480d0321e0SJeremy L Thompson 
490d0321e0SJeremy L Thompson   // Restrict
50437930d1SJeremy L Thompson   if (t_mode == CEED_NOTRANSPOSE) {
510d0321e0SJeremy L Thompson     // L-vector -> E-vector
520d0321e0SJeremy L Thompson     if (impl->d_ind) {
530d0321e0SJeremy L Thompson       // -- Offsets provided
54437930d1SJeremy L Thompson       kernel             = impl->OffsetNoTranspose;
55437930d1SJeremy L Thompson       void   *args[]     = {&num_elem, &impl->d_ind, &d_u, &d_v};
56437930d1SJeremy L Thompson       CeedInt block_size = elem_size < 256 ? (elem_size > 64 ? elem_size : 64) : 256;
572b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelHip(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
580d0321e0SJeremy L Thompson     } else {
590d0321e0SJeremy L Thompson       // -- Strided restriction
60437930d1SJeremy L Thompson       kernel             = impl->StridedNoTranspose;
61437930d1SJeremy L Thompson       void   *args[]     = {&num_elem, &d_u, &d_v};
62437930d1SJeremy L Thompson       CeedInt block_size = elem_size < 256 ? (elem_size > 64 ? elem_size : 64) : 256;
632b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelHip(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
640d0321e0SJeremy L Thompson     }
650d0321e0SJeremy L Thompson   } else {
660d0321e0SJeremy L Thompson     // E-vector -> L-vector
670d0321e0SJeremy L Thompson     if (impl->d_ind) {
680d0321e0SJeremy L Thompson       // -- Offsets provided
69437930d1SJeremy L Thompson       kernel       = impl->OffsetTranspose;
702b730f8bSJeremy L Thompson       void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
712b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelHip(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
720d0321e0SJeremy L Thompson     } else {
730d0321e0SJeremy L Thompson       // -- Strided restriction
74437930d1SJeremy L Thompson       kernel       = impl->StridedTranspose;
75437930d1SJeremy L Thompson       void *args[] = {&num_elem, &d_u, &d_v};
762b730f8bSJeremy L Thompson       CeedCallBackend(CeedRunKernelHip(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
770d0321e0SJeremy L Thompson     }
780d0321e0SJeremy L Thompson   }
790d0321e0SJeremy L Thompson 
802b730f8bSJeremy L Thompson   if (request != CEED_REQUEST_IMMEDIATE && request != CEED_REQUEST_ORDERED) *request = NULL;
810d0321e0SJeremy L Thompson 
820d0321e0SJeremy L Thompson   // Restore arrays
832b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
842b730f8bSJeremy L Thompson   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
850d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
860d0321e0SJeremy L Thompson }
870d0321e0SJeremy L Thompson 
880d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
890d0321e0SJeremy L Thompson // Blocked not supported
900d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
912b730f8bSJeremy L Thompson int CeedElemRestrictionApplyBlock_Hip(CeedElemRestriction r, CeedInt block, CeedTransposeMode t_mode, CeedVector u, CeedVector v,
922b730f8bSJeremy L Thompson                                       CeedRequest *request) {
930d0321e0SJeremy L Thompson   // LCOV_EXCL_START
940d0321e0SJeremy L Thompson   Ceed ceed;
952b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
962b730f8bSJeremy L Thompson   return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement blocked restrictions");
970d0321e0SJeremy L Thompson   // LCOV_EXCL_STOP
980d0321e0SJeremy L Thompson }
990d0321e0SJeremy L Thompson 
1000d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1010d0321e0SJeremy L Thompson // Get offsets
1020d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1032b730f8bSJeremy L Thompson static int CeedElemRestrictionGetOffsets_Hip(CeedElemRestriction rstr, CeedMemType mtype, const CeedInt **offsets) {
1040d0321e0SJeremy L Thompson   CeedElemRestriction_Hip *impl;
1052b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
1060d0321e0SJeremy L Thompson 
1070d0321e0SJeremy L Thompson   switch (mtype) {
1080d0321e0SJeremy L Thompson     case CEED_MEM_HOST:
1090d0321e0SJeremy L Thompson       *offsets = impl->h_ind;
1100d0321e0SJeremy L Thompson       break;
1110d0321e0SJeremy L Thompson     case CEED_MEM_DEVICE:
1120d0321e0SJeremy L Thompson       *offsets = impl->d_ind;
1130d0321e0SJeremy L Thompson       break;
1140d0321e0SJeremy L Thompson   }
1150d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1160d0321e0SJeremy L Thompson }
1170d0321e0SJeremy L Thompson 
1180d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1190d0321e0SJeremy L Thompson // Destroy restriction
1200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1210d0321e0SJeremy L Thompson static int CeedElemRestrictionDestroy_Hip(CeedElemRestriction r) {
1220d0321e0SJeremy L Thompson   CeedElemRestriction_Hip *impl;
1232b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(r, &impl));
1240d0321e0SJeremy L Thompson 
1250d0321e0SJeremy L Thompson   Ceed ceed;
1262b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
1272b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipModuleUnload(impl->module));
1282b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl->h_ind_allocated));
1292b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipFree(impl->d_ind_allocated));
1302b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipFree(impl->d_t_offsets));
1312b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipFree(impl->d_t_indices));
1322b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipFree(impl->d_l_vec_indices));
1332b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
134437930d1SJeremy L Thompson 
1350d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1360d0321e0SJeremy L Thompson }
1370d0321e0SJeremy L Thompson 
1380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1390d0321e0SJeremy L Thompson // Create transpose offsets and indices
1400d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1412b730f8bSJeremy L Thompson static int CeedElemRestrictionOffset_Hip(const CeedElemRestriction r, const CeedInt *indices) {
1420d0321e0SJeremy L Thompson   Ceed ceed;
1432b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
1440d0321e0SJeremy L Thompson   CeedElemRestriction_Hip *impl;
1452b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(r, &impl));
146e79b91d9SJeremy L Thompson   CeedSize l_size;
147e79b91d9SJeremy L Thompson   CeedInt  num_elem, elem_size, num_comp;
1482b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumElements(r, &num_elem));
1492b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size));
1502b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetLVectorSize(r, &l_size));
1512b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumComponents(r, &num_comp));
1520d0321e0SJeremy L Thompson 
153437930d1SJeremy L Thompson   // Count num_nodes
154437930d1SJeremy L Thompson   bool *is_node;
1552b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(l_size, &is_node));
156437930d1SJeremy L Thompson   const CeedInt size_indices = num_elem * elem_size;
1572b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < size_indices; i++) is_node[indices[i]] = 1;
158437930d1SJeremy L Thompson   CeedInt num_nodes = 0;
1592b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < l_size; i++) num_nodes += is_node[i];
160437930d1SJeremy L Thompson   impl->num_nodes = num_nodes;
1610d0321e0SJeremy L Thompson 
1620d0321e0SJeremy L Thompson   // L-vector offsets array
163437930d1SJeremy L Thompson   CeedInt *ind_to_offset, *l_vec_indices;
1642b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(l_size, &ind_to_offset));
1652b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(num_nodes, &l_vec_indices));
1660d0321e0SJeremy L Thompson   CeedInt j = 0;
1672b730f8bSJeremy L Thompson   for (CeedInt i = 0; i < l_size; i++) {
168437930d1SJeremy L Thompson     if (is_node[i]) {
169437930d1SJeremy L Thompson       l_vec_indices[j] = i;
1700d0321e0SJeremy L Thompson       ind_to_offset[i] = j++;
1710d0321e0SJeremy L Thompson     }
1722b730f8bSJeremy L Thompson   }
1732b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&is_node));
1740d0321e0SJeremy L Thompson 
1750d0321e0SJeremy L Thompson   // Compute transpose offsets and indices
176437930d1SJeremy L Thompson   const CeedInt size_offsets = num_nodes + 1;
177437930d1SJeremy L Thompson   CeedInt      *t_offsets;
1782b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(size_offsets, &t_offsets));
179437930d1SJeremy L Thompson   CeedInt *t_indices;
1802b730f8bSJeremy L Thompson   CeedCallBackend(CeedMalloc(size_indices, &t_indices));
1810d0321e0SJeremy L Thompson   // Count node multiplicity
1822b730f8bSJeremy L Thompson   for (CeedInt e = 0; e < num_elem; ++e) {
1832b730f8bSJeremy L Thompson     for (CeedInt i = 0; i < elem_size; ++i) ++t_offsets[ind_to_offset[indices[elem_size * e + i]] + 1];
1842b730f8bSJeremy L Thompson   }
1850d0321e0SJeremy L Thompson   // Convert to running sum
1862b730f8bSJeremy L Thompson   for (CeedInt i = 1; i < size_offsets; ++i) t_offsets[i] += t_offsets[i - 1];
1870d0321e0SJeremy L Thompson   // List all E-vec indices associated with L-vec node
188437930d1SJeremy L Thompson   for (CeedInt e = 0; e < num_elem; ++e) {
189437930d1SJeremy L Thompson     for (CeedInt i = 0; i < elem_size; ++i) {
190437930d1SJeremy L Thompson       const CeedInt lid                          = elem_size * e + i;
1910d0321e0SJeremy L Thompson       const CeedInt gid                          = indices[lid];
192437930d1SJeremy L Thompson       t_indices[t_offsets[ind_to_offset[gid]]++] = lid;
1930d0321e0SJeremy L Thompson     }
1940d0321e0SJeremy L Thompson   }
1950d0321e0SJeremy L Thompson   // Reset running sum
1962b730f8bSJeremy L Thompson   for (int i = size_offsets - 1; i > 0; --i) t_offsets[i] = t_offsets[i - 1];
197437930d1SJeremy L Thompson   t_offsets[0] = 0;
1980d0321e0SJeremy L Thompson 
1990d0321e0SJeremy L Thompson   // Copy data to device
2000d0321e0SJeremy L Thompson   // -- L-vector indices
2012b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&impl->d_l_vec_indices, num_nodes * sizeof(CeedInt)));
2022b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(impl->d_l_vec_indices, l_vec_indices, num_nodes * sizeof(CeedInt), hipMemcpyHostToDevice));
2030d0321e0SJeremy L Thompson   // -- Transpose offsets
2042b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&impl->d_t_offsets, size_offsets * sizeof(CeedInt)));
2052b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(impl->d_t_offsets, t_offsets, size_offsets * sizeof(CeedInt), hipMemcpyHostToDevice));
2060d0321e0SJeremy L Thompson   // -- Transpose indices
2072b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipMalloc((void **)&impl->d_t_indices, size_indices * sizeof(CeedInt)));
2082b730f8bSJeremy L Thompson   CeedCallHip(ceed, hipMemcpy(impl->d_t_indices, t_indices, size_indices * sizeof(CeedInt), hipMemcpyHostToDevice));
2090d0321e0SJeremy L Thompson 
2100d0321e0SJeremy L Thompson   // Cleanup
2112b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&ind_to_offset));
2122b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&l_vec_indices));
2132b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&t_offsets));
2142b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&t_indices));
215437930d1SJeremy L Thompson 
2160d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2170d0321e0SJeremy L Thompson }
2180d0321e0SJeremy L Thompson 
2190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2200d0321e0SJeremy L Thompson // Create restriction
2210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
2222b730f8bSJeremy L Thompson int CeedElemRestrictionCreate_Hip(CeedMemType mtype, CeedCopyMode cmode, const CeedInt *indices, CeedElemRestriction r) {
2230d0321e0SJeremy L Thompson   Ceed ceed;
2242b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
2250d0321e0SJeremy L Thompson   CeedElemRestriction_Hip *impl;
2262b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
227437930d1SJeremy L Thompson   CeedInt num_elem, num_comp, elem_size;
2282b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumElements(r, &num_elem));
2292b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumComponents(r, &num_comp));
2302b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size));
231437930d1SJeremy L Thompson   CeedInt size        = num_elem * elem_size;
232437930d1SJeremy L Thompson   CeedInt strides[3]  = {1, size, elem_size};
233437930d1SJeremy L Thompson   CeedInt comp_stride = 1;
2340d0321e0SJeremy L Thompson 
2350d0321e0SJeremy L Thompson   // Stride data
236437930d1SJeremy L Thompson   bool is_strided;
2372b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionIsStrided(r, &is_strided));
238437930d1SJeremy L Thompson   if (is_strided) {
239437930d1SJeremy L Thompson     bool has_backend_strides;
2402b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionHasBackendStrides(r, &has_backend_strides));
241437930d1SJeremy L Thompson     if (!has_backend_strides) {
2422b730f8bSJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetStrides(r, &strides));
2430d0321e0SJeremy L Thompson     }
2440d0321e0SJeremy L Thompson   } else {
2452b730f8bSJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(r, &comp_stride));
2460d0321e0SJeremy L Thompson   }
2470d0321e0SJeremy L Thompson 
2480d0321e0SJeremy L Thompson   impl->h_ind           = NULL;
2490d0321e0SJeremy L Thompson   impl->h_ind_allocated = NULL;
2500d0321e0SJeremy L Thompson   impl->d_ind           = NULL;
2510d0321e0SJeremy L Thompson   impl->d_ind_allocated = NULL;
252437930d1SJeremy L Thompson   impl->d_t_indices     = NULL;
253437930d1SJeremy L Thompson   impl->d_t_offsets     = NULL;
254437930d1SJeremy L Thompson   impl->num_nodes       = size;
2552b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionSetData(r, impl));
256437930d1SJeremy L Thompson   CeedInt layout[3] = {1, elem_size * num_elem, elem_size};
2572b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionSetELayout(r, layout));
2580d0321e0SJeremy L Thompson 
2590d0321e0SJeremy L Thompson   // Set up device indices/offset arrays
260*6574a04fSJeremy L Thompson   switch (mtype) {
261*6574a04fSJeremy L Thompson     case CEED_MEM_HOST: {
2620d0321e0SJeremy L Thompson       switch (cmode) {
2630d0321e0SJeremy L Thompson         case CEED_OWN_POINTER:
2640d0321e0SJeremy L Thompson           impl->h_ind_allocated = (CeedInt *)indices;
2650d0321e0SJeremy L Thompson           impl->h_ind           = (CeedInt *)indices;
2660d0321e0SJeremy L Thompson           break;
2670d0321e0SJeremy L Thompson         case CEED_USE_POINTER:
2680d0321e0SJeremy L Thompson           impl->h_ind = (CeedInt *)indices;
2690d0321e0SJeremy L Thompson           break;
2700d0321e0SJeremy L Thompson         case CEED_COPY_VALUES:
27144d7a66cSJeremy L Thompson           if (indices != NULL) {
2722b730f8bSJeremy L Thompson             CeedCallBackend(CeedMalloc(elem_size * num_elem, &impl->h_ind_allocated));
27344d7a66cSJeremy L Thompson             memcpy(impl->h_ind_allocated, indices, elem_size * num_elem * sizeof(CeedInt));
27444d7a66cSJeremy L Thompson             impl->h_ind = impl->h_ind_allocated;
27544d7a66cSJeremy L Thompson           }
2760d0321e0SJeremy L Thompson           break;
2770d0321e0SJeremy L Thompson       }
2780d0321e0SJeremy L Thompson       if (indices != NULL) {
2792b730f8bSJeremy L Thompson         CeedCallHip(ceed, hipMalloc((void **)&impl->d_ind, size * sizeof(CeedInt)));
2800d0321e0SJeremy L Thompson         impl->d_ind_allocated = impl->d_ind;  // We own the device memory
2812b730f8bSJeremy L Thompson         CeedCallHip(ceed, hipMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), hipMemcpyHostToDevice));
2822b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionOffset_Hip(r, indices));
2830d0321e0SJeremy L Thompson       }
284*6574a04fSJeremy L Thompson       break;
285*6574a04fSJeremy L Thompson     }
286*6574a04fSJeremy L Thompson     case CEED_MEM_DEVICE: {
2870d0321e0SJeremy L Thompson       switch (cmode) {
2880d0321e0SJeremy L Thompson         case CEED_COPY_VALUES:
2890d0321e0SJeremy L Thompson           if (indices != NULL) {
2902b730f8bSJeremy L Thompson             CeedCallHip(ceed, hipMalloc((void **)&impl->d_ind, size * sizeof(CeedInt)));
2910d0321e0SJeremy L Thompson             impl->d_ind_allocated = impl->d_ind;  // We own the device memory
2922b730f8bSJeremy L Thompson             CeedCallHip(ceed, hipMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), hipMemcpyDeviceToDevice));
2930d0321e0SJeremy L Thompson           }
2940d0321e0SJeremy L Thompson           break;
2950d0321e0SJeremy L Thompson         case CEED_OWN_POINTER:
2960d0321e0SJeremy L Thompson           impl->d_ind           = (CeedInt *)indices;
2970d0321e0SJeremy L Thompson           impl->d_ind_allocated = impl->d_ind;
2980d0321e0SJeremy L Thompson           break;
2990d0321e0SJeremy L Thompson         case CEED_USE_POINTER:
3000d0321e0SJeremy L Thompson           impl->d_ind = (CeedInt *)indices;
3010d0321e0SJeremy L Thompson       }
3020d0321e0SJeremy L Thompson       if (indices != NULL) {
3032b730f8bSJeremy L Thompson         CeedCallBackend(CeedMalloc(elem_size * num_elem, &impl->h_ind_allocated));
3042b730f8bSJeremy L Thompson         CeedCallHip(ceed, hipMemcpy(impl->h_ind_allocated, impl->d_ind, elem_size * num_elem * sizeof(CeedInt), hipMemcpyDeviceToHost));
30544d7a66cSJeremy L Thompson         impl->h_ind = impl->h_ind_allocated;
3062b730f8bSJeremy L Thompson         CeedCallBackend(CeedElemRestrictionOffset_Hip(r, indices));
3070d0321e0SJeremy L Thompson       }
308*6574a04fSJeremy L Thompson       break;
309*6574a04fSJeremy L Thompson     }
3100d0321e0SJeremy L Thompson     // LCOV_EXCL_START
311*6574a04fSJeremy L Thompson     default:
3122b730f8bSJeremy L Thompson       return CeedError(ceed, CEED_ERROR_BACKEND, "Only MemType = HOST or DEVICE supported");
3130d0321e0SJeremy L Thompson       // LCOV_EXCL_STOP
3140d0321e0SJeremy L Thompson   }
3150d0321e0SJeremy L Thompson 
3160d0321e0SJeremy L Thompson   // Compile HIP kernels
317437930d1SJeremy L Thompson   CeedInt num_nodes = impl->num_nodes;
318437930d1SJeremy L Thompson   char   *restriction_kernel_path, *restriction_kernel_source;
3192b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-restriction.h", &restriction_kernel_path));
32046dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Restriction Kernel Source -----\n");
3212b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source));
3222b730f8bSJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading Restriction Kernel Source Complete! -----\n");
3232b730f8bSJeremy L Thompson   CeedCallBackend(CeedCompileHip(ceed, restriction_kernel_source, &impl->module, 8, "RESTR_ELEM_SIZE", elem_size, "RESTR_NUM_ELEM", num_elem,
3242b730f8bSJeremy L Thompson                                  "RESTR_NUM_COMP", num_comp, "RESTR_NUM_NODES", num_nodes, "RESTR_COMP_STRIDE", comp_stride, "RESTR_STRIDE_NODES",
3252b730f8bSJeremy L Thompson                                  strides[0], "RESTR_STRIDE_COMP", strides[1], "RESTR_STRIDE_ELEM", strides[2]));
3262b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelHip(ceed, impl->module, "StridedNoTranspose", &impl->StridedNoTranspose));
3272b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelHip(ceed, impl->module, "OffsetNoTranspose", &impl->OffsetNoTranspose));
3282b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelHip(ceed, impl->module, "StridedTranspose", &impl->StridedTranspose));
3292b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelHip(ceed, impl->module, "OffsetTranspose", &impl->OffsetTranspose));
3302b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&restriction_kernel_path));
3312b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&restriction_kernel_source));
3320d0321e0SJeremy L Thompson 
3330d0321e0SJeremy L Thompson   // Register backend functions
3342b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "Apply", CeedElemRestrictionApply_Hip));
3352b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "ApplyBlock", CeedElemRestrictionApplyBlock_Hip));
3362b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "GetOffsets", CeedElemRestrictionGetOffsets_Hip));
3372b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "Destroy", CeedElemRestrictionDestroy_Hip));
3380d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3390d0321e0SJeremy L Thompson }
3400d0321e0SJeremy L Thompson 
3410d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3420d0321e0SJeremy L Thompson // Blocked not supported
3430d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
3442b730f8bSJeremy L Thompson int CeedElemRestrictionCreateBlocked_Hip(const CeedMemType mtype, const CeedCopyMode cmode, const CeedInt *indices, CeedElemRestriction r) {
3450d0321e0SJeremy L Thompson   Ceed ceed;
3462b730f8bSJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
3472b730f8bSJeremy L Thompson   return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement blocked restrictions");
3480d0321e0SJeremy L Thompson }
3490d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
350