xref: /libCEED/backends/hip-ref/ceed-hip-ref-restriction.c (revision 74d1b6c4d74ab01eefbd8940b6e3a2db84651b5c)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <stdbool.h>
12 #include <stddef.h>
13 #include <string.h>
14 #include <hip/hip_runtime.h>
15 
16 #include "../hip/ceed-hip-common.h"
17 #include "../hip/ceed-hip-compile.h"
18 #include "ceed-hip-ref.h"
19 
20 //------------------------------------------------------------------------------
21 // Core apply restriction code
22 //------------------------------------------------------------------------------
23 static inline int CeedElemRestrictionApply_Hip_Core(CeedElemRestriction rstr, CeedTransposeMode t_mode, bool use_signs, bool use_orients,
24                                                     CeedVector u, CeedVector v, CeedRequest *request) {
25   Ceed                     ceed;
26   CeedInt                  num_elem, elem_size;
27   CeedRestrictionType      rstr_type;
28   const CeedScalar        *d_u;
29   CeedScalar              *d_v;
30   CeedElemRestriction_Hip *impl;
31   hipFunction_t            kernel;
32 
33   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
34   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
35   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem));
36   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
37   CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type));
38   const CeedInt num_nodes = impl->num_nodes;
39 
40   // Get vectors
41   CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
42   if (t_mode == CEED_TRANSPOSE) {
43     // Sum into for transpose mode, e-vec to l-vec
44     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
45   } else {
46     // Overwrite for notranspose mode, l-vec to e-vec
47     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
48   }
49 
50   // Restrict
51   if (t_mode == CEED_NOTRANSPOSE) {
52     // L-vector -> E-vector
53     const CeedInt block_size = elem_size < 256 ? (elem_size > 64 ? elem_size : 64) : 256;
54     const CeedInt grid       = CeedDivUpInt(num_nodes, block_size);
55 
56     switch (rstr_type) {
57       case CEED_RESTRICTION_STRIDED: {
58         kernel       = impl->StridedNoTranspose;
59         void *args[] = {&num_elem, &d_u, &d_v};
60 
61         CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
62       } break;
63       case CEED_RESTRICTION_STANDARD: {
64         kernel       = impl->OffsetNoTranspose;
65         void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v};
66 
67         CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
68       } break;
69       case CEED_RESTRICTION_ORIENTED: {
70         if (use_signs) {
71           kernel       = impl->OrientedNoTranspose;
72           void *args[] = {&num_elem, &impl->d_ind, &impl->d_orients, &d_u, &d_v};
73 
74           CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
75         } else {
76           kernel       = impl->OffsetNoTranspose;
77           void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v};
78 
79           CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
80         }
81       } break;
82       case CEED_RESTRICTION_CURL_ORIENTED: {
83         if (use_signs && use_orients) {
84           kernel       = impl->CurlOrientedNoTranspose;
85           void *args[] = {&num_elem, &impl->d_ind, &impl->d_curl_orients, &d_u, &d_v};
86 
87           CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
88         } else if (use_orients) {
89           kernel       = impl->CurlOrientedUnsignedNoTranspose;
90           void *args[] = {&num_elem, &impl->d_ind, &impl->d_curl_orients, &d_u, &d_v};
91 
92           CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
93         } else {
94           kernel       = impl->OffsetNoTranspose;
95           void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v};
96 
97           CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
98         }
99       } break;
100       case CEED_RESTRICTION_POINTS: {
101         // LCOV_EXCL_START
102         return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement restriction CeedElemRestrictionAtPoints");
103         // LCOV_EXCL_STOP
104       } break;
105     }
106   } else {
107     // E-vector -> L-vector
108     const CeedInt block_size = 64;
109     const CeedInt grid       = CeedDivUpInt(num_nodes, block_size);
110 
111     switch (rstr_type) {
112       case CEED_RESTRICTION_STRIDED: {
113         kernel       = impl->StridedTranspose;
114         void *args[] = {&num_elem, &d_u, &d_v};
115 
116         CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
117       } break;
118       case CEED_RESTRICTION_STANDARD: {
119         if (impl->OffsetTranspose) {
120           kernel       = impl->OffsetTranspose;
121           void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v};
122 
123           CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
124         } else {
125           kernel       = impl->OffsetTransposeDet;
126           void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
127 
128           CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
129         }
130       } break;
131       case CEED_RESTRICTION_ORIENTED: {
132         if (use_signs) {
133           if (impl->OrientedTranspose) {
134             kernel       = impl->OrientedTranspose;
135             void *args[] = {&num_elem, &impl->d_ind, &impl->d_orients, &d_u, &d_v};
136 
137             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
138           } else {
139             kernel       = impl->OrientedTransposeDet;
140             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_orients, &d_u, &d_v};
141 
142             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
143           }
144         } else {
145           if (impl->OffsetTranspose) {
146             kernel       = impl->OffsetTranspose;
147             void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v};
148 
149             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
150           } else {
151             kernel       = impl->OffsetTransposeDet;
152             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
153 
154             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
155           }
156         }
157       } break;
158       case CEED_RESTRICTION_CURL_ORIENTED: {
159         if (use_signs && use_orients) {
160           if (impl->CurlOrientedTranspose) {
161             kernel       = impl->CurlOrientedTranspose;
162             void *args[] = {&num_elem, &impl->d_ind, &impl->d_curl_orients, &d_u, &d_v};
163 
164             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
165           } else {
166             kernel       = impl->CurlOrientedTransposeDet;
167             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v};
168 
169             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
170           }
171         } else if (use_orients) {
172           if (impl->CurlOrientedUnsignedTranspose) {
173             kernel       = impl->CurlOrientedUnsignedTranspose;
174             void *args[] = {&num_elem, &impl->d_ind, &impl->d_curl_orients, &d_u, &d_v};
175 
176             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
177           } else {
178             kernel       = impl->CurlOrientedUnsignedTransposeDet;
179             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v};
180 
181             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
182           }
183         } else {
184           if (impl->OffsetTranspose) {
185             kernel       = impl->OffsetTranspose;
186             void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v};
187 
188             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
189           } else {
190             kernel       = impl->OffsetTransposeDet;
191             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
192 
193             CeedCallBackend(CeedRunKernel_Hip(ceed, kernel, grid, block_size, args));
194           }
195         }
196       } break;
197       case CEED_RESTRICTION_POINTS: {
198         // LCOV_EXCL_START
199         return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement restriction CeedElemRestrictionAtPoints");
200         // LCOV_EXCL_STOP
201       } break;
202     }
203   }
204 
205   if (request != CEED_REQUEST_IMMEDIATE && request != CEED_REQUEST_ORDERED) *request = NULL;
206 
207   // Restore arrays
208   CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
209   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
210   return CEED_ERROR_SUCCESS;
211 }
212 
213 //------------------------------------------------------------------------------
214 // Apply restriction
215 //------------------------------------------------------------------------------
216 static int CeedElemRestrictionApply_Hip(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, CeedRequest *request) {
217   return CeedElemRestrictionApply_Hip_Core(rstr, t_mode, true, true, u, v, request);
218 }
219 
220 //------------------------------------------------------------------------------
221 // Apply unsigned restriction
222 //------------------------------------------------------------------------------
223 static int CeedElemRestrictionApplyUnsigned_Hip(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v,
224                                                 CeedRequest *request) {
225   return CeedElemRestrictionApply_Hip_Core(rstr, t_mode, false, true, u, v, request);
226 }
227 
228 //------------------------------------------------------------------------------
229 // Apply unoriented restriction
230 //------------------------------------------------------------------------------
231 static int CeedElemRestrictionApplyUnoriented_Hip(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v,
232                                                   CeedRequest *request) {
233   return CeedElemRestrictionApply_Hip_Core(rstr, t_mode, false, false, u, v, request);
234 }
235 
236 //------------------------------------------------------------------------------
237 // Get offsets
238 //------------------------------------------------------------------------------
239 static int CeedElemRestrictionGetOffsets_Hip(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt **offsets) {
240   CeedElemRestriction_Hip *impl;
241 
242   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
243   switch (mem_type) {
244     case CEED_MEM_HOST:
245       *offsets = impl->h_ind;
246       break;
247     case CEED_MEM_DEVICE:
248       *offsets = impl->d_ind;
249       break;
250   }
251   return CEED_ERROR_SUCCESS;
252 }
253 
254 //------------------------------------------------------------------------------
255 // Get orientations
256 //------------------------------------------------------------------------------
257 static int CeedElemRestrictionGetOrientations_Hip(CeedElemRestriction rstr, CeedMemType mem_type, const bool **orients) {
258   CeedElemRestriction_Hip *impl;
259   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
260 
261   switch (mem_type) {
262     case CEED_MEM_HOST:
263       *orients = impl->h_orients;
264       break;
265     case CEED_MEM_DEVICE:
266       *orients = impl->d_orients;
267       break;
268   }
269   return CEED_ERROR_SUCCESS;
270 }
271 
272 //------------------------------------------------------------------------------
273 // Get curl-conforming orientations
274 //------------------------------------------------------------------------------
275 static int CeedElemRestrictionGetCurlOrientations_Hip(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt8 **curl_orients) {
276   CeedElemRestriction_Hip *impl;
277   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
278 
279   switch (mem_type) {
280     case CEED_MEM_HOST:
281       *curl_orients = impl->h_curl_orients;
282       break;
283     case CEED_MEM_DEVICE:
284       *curl_orients = impl->d_curl_orients;
285       break;
286   }
287   return CEED_ERROR_SUCCESS;
288 }
289 
290 //------------------------------------------------------------------------------
291 // Destroy restriction
292 //------------------------------------------------------------------------------
293 static int CeedElemRestrictionDestroy_Hip(CeedElemRestriction rstr) {
294   Ceed                     ceed;
295   CeedElemRestriction_Hip *impl;
296 
297   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
298   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
299   CeedCallHip(ceed, hipModuleUnload(impl->module));
300   CeedCallBackend(CeedFree(&impl->h_ind_allocated));
301   CeedCallHip(ceed, hipFree(impl->d_ind_allocated));
302   CeedCallHip(ceed, hipFree(impl->d_t_offsets));
303   CeedCallHip(ceed, hipFree(impl->d_t_indices));
304   CeedCallHip(ceed, hipFree(impl->d_l_vec_indices));
305   CeedCallBackend(CeedFree(&impl->h_orients_allocated));
306   CeedCallHip(ceed, hipFree(impl->d_orients_allocated));
307   CeedCallBackend(CeedFree(&impl->h_curl_orients_allocated));
308   CeedCallHip(ceed, hipFree(impl->d_curl_orients_allocated));
309   CeedCallBackend(CeedFree(&impl));
310   return CEED_ERROR_SUCCESS;
311 }
312 
313 //------------------------------------------------------------------------------
314 // Create transpose offsets and indices
315 //------------------------------------------------------------------------------
316 static int CeedElemRestrictionOffset_Hip(const CeedElemRestriction rstr, const CeedInt *indices) {
317   Ceed                     ceed;
318   bool                    *is_node;
319   CeedSize                 l_size;
320   CeedInt                  num_elem, elem_size, num_comp, num_nodes = 0;
321   CeedInt                 *ind_to_offset, *l_vec_indices, *t_offsets, *t_indices;
322   CeedElemRestriction_Hip *impl;
323 
324   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
325   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
326   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem));
327   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
328   CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size));
329   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp));
330   const CeedInt size_indices = num_elem * elem_size;
331 
332   // Count num_nodes
333   CeedCallBackend(CeedCalloc(l_size, &is_node));
334 
335   for (CeedInt i = 0; i < size_indices; i++) is_node[indices[i]] = 1;
336   for (CeedInt i = 0; i < l_size; i++) num_nodes += is_node[i];
337   impl->num_nodes = num_nodes;
338 
339   // L-vector offsets array
340   CeedCallBackend(CeedCalloc(l_size, &ind_to_offset));
341   CeedCallBackend(CeedCalloc(num_nodes, &l_vec_indices));
342   for (CeedInt i = 0, j = 0; i < l_size; i++) {
343     if (is_node[i]) {
344       l_vec_indices[j] = i;
345       ind_to_offset[i] = j++;
346     }
347   }
348   CeedCallBackend(CeedFree(&is_node));
349 
350   // Compute transpose offsets and indices
351   const CeedInt size_offsets = num_nodes + 1;
352 
353   CeedCallBackend(CeedCalloc(size_offsets, &t_offsets));
354   CeedCallBackend(CeedMalloc(size_indices, &t_indices));
355   // Count node multiplicity
356   for (CeedInt e = 0; e < num_elem; ++e) {
357     for (CeedInt i = 0; i < elem_size; ++i) ++t_offsets[ind_to_offset[indices[elem_size * e + i]] + 1];
358   }
359   // Convert to running sum
360   for (CeedInt i = 1; i < size_offsets; ++i) t_offsets[i] += t_offsets[i - 1];
361   // List all E-vec indices associated with L-vec node
362   for (CeedInt e = 0; e < num_elem; ++e) {
363     for (CeedInt i = 0; i < elem_size; ++i) {
364       const CeedInt lid = elem_size * e + i;
365       const CeedInt gid = indices[lid];
366 
367       t_indices[t_offsets[ind_to_offset[gid]]++] = lid;
368     }
369   }
370   // Reset running sum
371   for (int i = size_offsets - 1; i > 0; --i) t_offsets[i] = t_offsets[i - 1];
372   t_offsets[0] = 0;
373 
374   // Copy data to device
375   // -- L-vector indices
376   CeedCallHip(ceed, hipMalloc((void **)&impl->d_l_vec_indices, num_nodes * sizeof(CeedInt)));
377   CeedCallHip(ceed, hipMemcpy(impl->d_l_vec_indices, l_vec_indices, num_nodes * sizeof(CeedInt), hipMemcpyHostToDevice));
378   // -- Transpose offsets
379   CeedCallHip(ceed, hipMalloc((void **)&impl->d_t_offsets, size_offsets * sizeof(CeedInt)));
380   CeedCallHip(ceed, hipMemcpy(impl->d_t_offsets, t_offsets, size_offsets * sizeof(CeedInt), hipMemcpyHostToDevice));
381   // -- Transpose indices
382   CeedCallHip(ceed, hipMalloc((void **)&impl->d_t_indices, size_indices * sizeof(CeedInt)));
383   CeedCallHip(ceed, hipMemcpy(impl->d_t_indices, t_indices, size_indices * sizeof(CeedInt), hipMemcpyHostToDevice));
384 
385   // Cleanup
386   CeedCallBackend(CeedFree(&ind_to_offset));
387   CeedCallBackend(CeedFree(&l_vec_indices));
388   CeedCallBackend(CeedFree(&t_offsets));
389   CeedCallBackend(CeedFree(&t_indices));
390   return CEED_ERROR_SUCCESS;
391 }
392 
393 //------------------------------------------------------------------------------
394 // Create restriction
395 //------------------------------------------------------------------------------
396 int CeedElemRestrictionCreate_Hip(CeedMemType mem_type, CeedCopyMode copy_mode, const CeedInt *indices, const bool *orients,
397                                   const CeedInt8 *curl_orients, CeedElemRestriction rstr) {
398   Ceed                     ceed, ceed_parent;
399   bool                     is_deterministic;
400   CeedInt                  num_elem, num_comp, elem_size, comp_stride = 1;
401   CeedRestrictionType      rstr_type;
402   char                    *restriction_kernel_path, *restriction_kernel_source;
403   CeedElemRestriction_Hip *impl;
404 
405   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
406   CeedCallBackend(CeedGetParent(ceed, &ceed_parent));
407   CeedCallBackend(CeedIsDeterministic(ceed_parent, &is_deterministic));
408   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem));
409   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp));
410   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
411   const CeedInt size       = num_elem * elem_size;
412   CeedInt       strides[3] = {1, size, elem_size};
413   CeedInt       layout[3]  = {1, elem_size * num_elem, elem_size};
414 
415   // Stride data
416   CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type));
417   if (rstr_type == CEED_RESTRICTION_STRIDED) {
418     bool has_backend_strides;
419 
420     CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &has_backend_strides));
421     if (!has_backend_strides) {
422       CeedCallBackend(CeedElemRestrictionGetStrides(rstr, &strides));
423     }
424   } else {
425     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &comp_stride));
426   }
427 
428   CeedCallBackend(CeedCalloc(1, &impl));
429   impl->num_nodes                = size;
430   impl->h_ind                    = NULL;
431   impl->h_ind_allocated          = NULL;
432   impl->d_ind                    = NULL;
433   impl->d_ind_allocated          = NULL;
434   impl->d_t_indices              = NULL;
435   impl->d_t_offsets              = NULL;
436   impl->h_orients                = NULL;
437   impl->h_orients_allocated      = NULL;
438   impl->d_orients                = NULL;
439   impl->d_orients_allocated      = NULL;
440   impl->h_curl_orients           = NULL;
441   impl->h_curl_orients_allocated = NULL;
442   impl->d_curl_orients           = NULL;
443   impl->d_curl_orients_allocated = NULL;
444   CeedCallBackend(CeedElemRestrictionSetData(rstr, impl));
445   CeedCallBackend(CeedElemRestrictionSetELayout(rstr, layout));
446 
447   // Set up device offset/orientation arrays
448   if (rstr_type != CEED_RESTRICTION_STRIDED) {
449     switch (mem_type) {
450       case CEED_MEM_HOST: {
451         switch (copy_mode) {
452           case CEED_OWN_POINTER:
453             impl->h_ind_allocated = (CeedInt *)indices;
454             impl->h_ind           = (CeedInt *)indices;
455             break;
456           case CEED_USE_POINTER:
457             impl->h_ind = (CeedInt *)indices;
458             break;
459           case CEED_COPY_VALUES:
460             CeedCallBackend(CeedMalloc(size, &impl->h_ind_allocated));
461             memcpy(impl->h_ind_allocated, indices, size * sizeof(CeedInt));
462             impl->h_ind = impl->h_ind_allocated;
463             break;
464         }
465         CeedCallHip(ceed, hipMalloc((void **)&impl->d_ind, size * sizeof(CeedInt)));
466         impl->d_ind_allocated = impl->d_ind;  // We own the device memory
467         CeedCallHip(ceed, hipMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), hipMemcpyHostToDevice));
468         if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Hip(rstr, indices));
469       } break;
470       case CEED_MEM_DEVICE: {
471         switch (copy_mode) {
472           case CEED_COPY_VALUES:
473             CeedCallHip(ceed, hipMalloc((void **)&impl->d_ind, size * sizeof(CeedInt)));
474             impl->d_ind_allocated = impl->d_ind;  // We own the device memory
475             CeedCallHip(ceed, hipMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), hipMemcpyDeviceToDevice));
476             break;
477           case CEED_OWN_POINTER:
478             impl->d_ind           = (CeedInt *)indices;
479             impl->d_ind_allocated = impl->d_ind;
480             break;
481           case CEED_USE_POINTER:
482             impl->d_ind = (CeedInt *)indices;
483             break;
484         }
485         CeedCallBackend(CeedMalloc(size, &impl->h_ind_allocated));
486         CeedCallHip(ceed, hipMemcpy(impl->h_ind_allocated, impl->d_ind, size * sizeof(CeedInt), hipMemcpyDeviceToHost));
487         impl->h_ind = impl->h_ind_allocated;
488         if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Hip(rstr, indices));
489       } break;
490     }
491 
492     // Orientation data
493     if (rstr_type == CEED_RESTRICTION_ORIENTED) {
494       switch (mem_type) {
495         case CEED_MEM_HOST: {
496           switch (copy_mode) {
497             case CEED_OWN_POINTER:
498               impl->h_orients_allocated = (bool *)orients;
499               impl->h_orients           = (bool *)orients;
500               break;
501             case CEED_USE_POINTER:
502               impl->h_orients = (bool *)orients;
503               break;
504             case CEED_COPY_VALUES:
505               CeedCallBackend(CeedMalloc(size, &impl->h_orients_allocated));
506               memcpy(impl->h_orients_allocated, orients, size * sizeof(bool));
507               impl->h_orients = impl->h_orients_allocated;
508               break;
509           }
510           CeedCallHip(ceed, hipMalloc((void **)&impl->d_orients, size * sizeof(bool)));
511           impl->d_orients_allocated = impl->d_orients;  // We own the device memory
512           CeedCallHip(ceed, hipMemcpy(impl->d_orients, orients, size * sizeof(bool), hipMemcpyHostToDevice));
513         } break;
514         case CEED_MEM_DEVICE: {
515           switch (copy_mode) {
516             case CEED_COPY_VALUES:
517               CeedCallHip(ceed, hipMalloc((void **)&impl->d_orients, size * sizeof(bool)));
518               impl->d_orients_allocated = impl->d_orients;  // We own the device memory
519               CeedCallHip(ceed, hipMemcpy(impl->d_orients, orients, size * sizeof(bool), hipMemcpyDeviceToDevice));
520               break;
521             case CEED_OWN_POINTER:
522               impl->d_orients           = (bool *)orients;
523               impl->d_orients_allocated = impl->d_orients;
524               break;
525             case CEED_USE_POINTER:
526               impl->d_orients = (bool *)orients;
527               break;
528           }
529           CeedCallBackend(CeedMalloc(size, &impl->h_orients_allocated));
530           CeedCallHip(ceed, hipMemcpy(impl->h_orients_allocated, impl->d_orients, size * sizeof(bool), hipMemcpyDeviceToHost));
531           impl->h_orients = impl->h_orients_allocated;
532         } break;
533       }
534     } else if (rstr_type == CEED_RESTRICTION_CURL_ORIENTED) {
535       switch (mem_type) {
536         case CEED_MEM_HOST: {
537           switch (copy_mode) {
538             case CEED_OWN_POINTER:
539               impl->h_curl_orients_allocated = (CeedInt8 *)curl_orients;
540               impl->h_curl_orients           = (CeedInt8 *)curl_orients;
541               break;
542             case CEED_USE_POINTER:
543               impl->h_curl_orients = (CeedInt8 *)curl_orients;
544               break;
545             case CEED_COPY_VALUES:
546               CeedCallBackend(CeedMalloc(3 * size, &impl->h_curl_orients_allocated));
547               memcpy(impl->h_curl_orients_allocated, curl_orients, 3 * size * sizeof(CeedInt8));
548               impl->h_curl_orients = impl->h_curl_orients_allocated;
549               break;
550           }
551           CeedCallHip(ceed, hipMalloc((void **)&impl->d_curl_orients, 3 * size * sizeof(CeedInt8)));
552           impl->d_curl_orients_allocated = impl->d_curl_orients;  // We own the device memory
553           CeedCallHip(ceed, hipMemcpy(impl->d_curl_orients, curl_orients, 3 * size * sizeof(CeedInt8), hipMemcpyHostToDevice));
554         } break;
555         case CEED_MEM_DEVICE: {
556           switch (copy_mode) {
557             case CEED_COPY_VALUES:
558               CeedCallHip(ceed, hipMalloc((void **)&impl->d_curl_orients, 3 * size * sizeof(CeedInt8)));
559               impl->d_curl_orients_allocated = impl->d_curl_orients;  // We own the device memory
560               CeedCallHip(ceed, hipMemcpy(impl->d_curl_orients, curl_orients, 3 * size * sizeof(CeedInt8), hipMemcpyDeviceToDevice));
561               break;
562             case CEED_OWN_POINTER:
563               impl->d_curl_orients           = (CeedInt8 *)curl_orients;
564               impl->d_curl_orients_allocated = impl->d_curl_orients;
565               break;
566             case CEED_USE_POINTER:
567               impl->d_curl_orients = (CeedInt8 *)curl_orients;
568               break;
569           }
570           CeedCallBackend(CeedMalloc(3 * size, &impl->h_curl_orients_allocated));
571           CeedCallHip(ceed, hipMemcpy(impl->h_curl_orients_allocated, impl->d_curl_orients, 3 * size * sizeof(CeedInt8), hipMemcpyDeviceToHost));
572           impl->h_curl_orients = impl->h_curl_orients_allocated;
573         } break;
574       }
575     }
576   }
577 
578   // Compile HIP kernels
579   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-restriction.h", &restriction_kernel_path));
580   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n");
581   CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source));
582   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n");
583   CeedCallBackend(CeedCompile_Hip(ceed, restriction_kernel_source, &impl->module, 8, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem,
584                                   "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride, "RSTR_STRIDE_NODES",
585                                   strides[0], "RSTR_STRIDE_COMP", strides[1], "RSTR_STRIDE_ELEM", strides[2]));
586   CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "StridedNoTranspose", &impl->StridedNoTranspose));
587   CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "StridedTranspose", &impl->StridedTranspose));
588   CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "OffsetNoTranspose", &impl->OffsetNoTranspose));
589   CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "OrientedNoTranspose", &impl->OrientedNoTranspose));
590   CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "CurlOrientedNoTranspose", &impl->CurlOrientedNoTranspose));
591   CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "CurlOrientedUnsignedNoTranspose", &impl->CurlOrientedUnsignedNoTranspose));
592   if (!is_deterministic) {
593     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "OffsetTranspose", &impl->OffsetTranspose));
594     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "OrientedTranspose", &impl->OrientedTranspose));
595     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "CurlOrientedTranspose", &impl->CurlOrientedTranspose));
596     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "CurlOrientedUnsignedTranspose", &impl->CurlOrientedUnsignedTranspose));
597   } else {
598     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "OffsetTransposeDet", &impl->OffsetTransposeDet));
599     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "OrientedTransposeDet", &impl->OrientedTransposeDet));
600     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "CurlOrientedTransposeDet", &impl->CurlOrientedTransposeDet));
601     CeedCallBackend(CeedGetKernel_Hip(ceed, impl->module, "CurlOrientedUnsignedTransposeDet", &impl->CurlOrientedUnsignedTransposeDet));
602   }
603   CeedCallBackend(CeedFree(&restriction_kernel_path));
604   CeedCallBackend(CeedFree(&restriction_kernel_source));
605 
606   // Register backend functions
607   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Apply", CeedElemRestrictionApply_Hip));
608   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnsigned", CeedElemRestrictionApplyUnsigned_Hip));
609   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnoriented", CeedElemRestrictionApplyUnoriented_Hip));
610   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOffsets", CeedElemRestrictionGetOffsets_Hip));
611   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOrientations", CeedElemRestrictionGetOrientations_Hip));
612   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetCurlOrientations", CeedElemRestrictionGetCurlOrientations_Hip));
613   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Destroy", CeedElemRestrictionDestroy_Hip));
614   return CEED_ERROR_SUCCESS;
615 }
616 
617 //------------------------------------------------------------------------------
618