xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-restriction.c (revision 9c25dd66b9687765a7022cc762ccaf201b721845) !
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2ff1e7120SSebastian Grimberg // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ff1e7120SSebastian Grimberg //
4ff1e7120SSebastian Grimberg // SPDX-License-Identifier: BSD-2-Clause
5ff1e7120SSebastian Grimberg //
6ff1e7120SSebastian Grimberg // This file is part of CEED:  http://github.com/ceed
7ff1e7120SSebastian Grimberg 
8ff1e7120SSebastian Grimberg #include <ceed.h>
9ff1e7120SSebastian Grimberg #include <ceed/backend.h>
10ff1e7120SSebastian Grimberg #include <ceed/jit-tools.h>
11ff1e7120SSebastian Grimberg #include <cuda.h>
12ff1e7120SSebastian Grimberg #include <cuda_runtime.h>
13ff1e7120SSebastian Grimberg #include <stdbool.h>
14ff1e7120SSebastian Grimberg #include <stddef.h>
15ff1e7120SSebastian Grimberg #include <string.h>
16ff1e7120SSebastian Grimberg 
17ff1e7120SSebastian Grimberg #include "../cuda/ceed-cuda-common.h"
18ff1e7120SSebastian Grimberg #include "../cuda/ceed-cuda-compile.h"
19ff1e7120SSebastian Grimberg #include "ceed-cuda-ref.h"
20ff1e7120SSebastian Grimberg 
21ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
22cf8cbdd6SSebastian Grimberg // Compile restriction kernels
23cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
24cf8cbdd6SSebastian Grimberg static inline int CeedElemRestrictionSetupCompile_Cuda(CeedElemRestriction rstr) {
25cf8cbdd6SSebastian Grimberg   Ceed                      ceed;
26cf8cbdd6SSebastian Grimberg   bool                      is_deterministic;
27cf8cbdd6SSebastian Grimberg   CeedInt                   num_elem, num_comp, elem_size, comp_stride;
28cf8cbdd6SSebastian Grimberg   CeedRestrictionType       rstr_type;
29cf8cbdd6SSebastian Grimberg   CeedElemRestriction_Cuda *impl;
30cf8cbdd6SSebastian Grimberg 
31cf8cbdd6SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
32cf8cbdd6SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
33b20a4af9SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type));
34cf8cbdd6SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem));
35cf8cbdd6SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp));
36cf8cbdd6SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &comp_stride));
37b20a4af9SJeremy L Thompson   if (rstr_type == CEED_RESTRICTION_POINTS) {
38b20a4af9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr, &elem_size));
39b20a4af9SJeremy L Thompson   } else {
40b20a4af9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
41b20a4af9SJeremy L Thompson   }
42cf8cbdd6SSebastian Grimberg   is_deterministic = impl->d_l_vec_indices != NULL;
43cf8cbdd6SSebastian Grimberg 
44cf8cbdd6SSebastian Grimberg   // Compile CUDA kernels
45cf8cbdd6SSebastian Grimberg   switch (rstr_type) {
46cf8cbdd6SSebastian Grimberg     case CEED_RESTRICTION_STRIDED: {
47*9c25dd66SJeremy L Thompson       const char restriction_kernel_source[] = "// Strided restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-strided.h>\n";
48cf8cbdd6SSebastian Grimberg       bool       has_backend_strides;
49509d4af6SJeremy L Thompson       CeedInt    strides[3] = {1, num_elem * elem_size, elem_size};
50cf8cbdd6SSebastian Grimberg 
51cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &has_backend_strides));
52cf8cbdd6SSebastian Grimberg       if (!has_backend_strides) {
5356c48462SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(rstr, strides));
54cf8cbdd6SSebastian Grimberg       }
55cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem,
56cf8cbdd6SSebastian Grimberg                                        "RSTR_NUM_COMP", num_comp, "RSTR_STRIDE_NODES", strides[0], "RSTR_STRIDE_COMP", strides[1], "RSTR_STRIDE_ELEM",
57cf8cbdd6SSebastian Grimberg                                        strides[2]));
58cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedNoTranspose", &impl->ApplyNoTranspose));
59cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedTranspose", &impl->ApplyTranspose));
60cf8cbdd6SSebastian Grimberg     } break;
61cf8cbdd6SSebastian Grimberg     case CEED_RESTRICTION_STANDARD: {
62*9c25dd66SJeremy L Thompson       const char restriction_kernel_source[] = "// Standard restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-offset.h>\n";
63*9c25dd66SJeremy L Thompson 
64cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem,
65cf8cbdd6SSebastian Grimberg                                        "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride,
66cf8cbdd6SSebastian Grimberg                                        "USE_DETERMINISTIC", is_deterministic ? 1 : 0));
67cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyNoTranspose));
68cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyTranspose));
69cf8cbdd6SSebastian Grimberg     } break;
70*9c25dd66SJeremy L Thompson     case CEED_RESTRICTION_POINTS: {
71*9c25dd66SJeremy L Thompson       const char restriction_kernel_source[] =
72*9c25dd66SJeremy L Thompson           "// AtPoints restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-at-points.h>\n\n"
73*9c25dd66SJeremy L Thompson           "// Standard restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-offset.h>\n";
74cf8cbdd6SSebastian Grimberg 
75*9c25dd66SJeremy L Thompson       CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem,
76*9c25dd66SJeremy L Thompson                                        "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride,
77*9c25dd66SJeremy L Thompson                                        "USE_DETERMINISTIC", is_deterministic ? 1 : 0));
78*9c25dd66SJeremy L Thompson       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyNoTranspose));
79*9c25dd66SJeremy L Thompson       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "AtPointsTranspose", &impl->ApplyTranspose));
80*9c25dd66SJeremy L Thompson     } break;
81*9c25dd66SJeremy L Thompson     case CEED_RESTRICTION_ORIENTED: {
82*9c25dd66SJeremy L Thompson       const char restriction_kernel_source[] =
83*9c25dd66SJeremy L Thompson           "// Oriented restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-oriented.h>\n\n"
84*9c25dd66SJeremy L Thompson           "// Standard restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-offset.h>\n";
85*9c25dd66SJeremy L Thompson 
86cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem,
87cf8cbdd6SSebastian Grimberg                                        "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride,
88cf8cbdd6SSebastian Grimberg                                        "USE_DETERMINISTIC", is_deterministic ? 1 : 0));
89cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OrientedNoTranspose", &impl->ApplyNoTranspose));
90cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyUnsignedNoTranspose));
91cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OrientedTranspose", &impl->ApplyTranspose));
92cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyUnsignedTranspose));
93cf8cbdd6SSebastian Grimberg     } break;
94cf8cbdd6SSebastian Grimberg     case CEED_RESTRICTION_CURL_ORIENTED: {
95*9c25dd66SJeremy L Thompson       const char restriction_kernel_source[] =
96*9c25dd66SJeremy L Thompson           "// Curl oriented restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-curl-oriented.h>\n\n"
97*9c25dd66SJeremy L Thompson           "// Standard restriction source\n#include <ceed/jit-source/cuda/cuda-ref-restriction-offset.h>\n";
98cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem,
99cf8cbdd6SSebastian Grimberg                                        "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride,
100cf8cbdd6SSebastian Grimberg                                        "USE_DETERMINISTIC", is_deterministic ? 1 : 0));
101cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedNoTranspose", &impl->ApplyNoTranspose));
102cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedUnsignedNoTranspose", &impl->ApplyUnsignedNoTranspose));
103cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyUnorientedNoTranspose));
104cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedTranspose", &impl->ApplyTranspose));
105cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedUnsignedTranspose", &impl->ApplyUnsignedTranspose));
106cf8cbdd6SSebastian Grimberg       CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyUnorientedTranspose));
107cf8cbdd6SSebastian Grimberg     } break;
108cf8cbdd6SSebastian Grimberg   }
109cf8cbdd6SSebastian Grimberg   return CEED_ERROR_SUCCESS;
110cf8cbdd6SSebastian Grimberg }
111cf8cbdd6SSebastian Grimberg 
112cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------
113dce49693SSebastian Grimberg // Core apply restriction code
114ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
115dce49693SSebastian Grimberg static inline int CeedElemRestrictionApply_Cuda_Core(CeedElemRestriction rstr, CeedTransposeMode t_mode, bool use_signs, bool use_orients,
116dce49693SSebastian Grimberg                                                      CeedVector u, CeedVector v, CeedRequest *request) {
117ff1e7120SSebastian Grimberg   Ceed                      ceed;
118dce49693SSebastian Grimberg   CeedRestrictionType       rstr_type;
119ff1e7120SSebastian Grimberg   const CeedScalar         *d_u;
120ff1e7120SSebastian Grimberg   CeedScalar               *d_v;
121ca735530SJeremy L Thompson   CeedElemRestriction_Cuda *impl;
122ca735530SJeremy L Thompson 
123dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
124dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
125dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type));
126cf8cbdd6SSebastian Grimberg 
127cf8cbdd6SSebastian Grimberg   // Assemble kernel if needed
128cf8cbdd6SSebastian Grimberg   if (!impl->module) {
129cf8cbdd6SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionSetupCompile_Cuda(rstr));
130cf8cbdd6SSebastian Grimberg   }
131ca735530SJeremy L Thompson 
132ca735530SJeremy L Thompson   // Get vectors
133ff1e7120SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
134ff1e7120SSebastian Grimberg   if (t_mode == CEED_TRANSPOSE) {
135ff1e7120SSebastian Grimberg     // Sum into for transpose mode, e-vec to l-vec
136ff1e7120SSebastian Grimberg     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
137ff1e7120SSebastian Grimberg   } else {
138ff1e7120SSebastian Grimberg     // Overwrite for notranspose mode, l-vec to e-vec
139ff1e7120SSebastian Grimberg     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
140ff1e7120SSebastian Grimberg   }
141ff1e7120SSebastian Grimberg 
142ff1e7120SSebastian Grimberg   // Restrict
143ff1e7120SSebastian Grimberg   if (t_mode == CEED_NOTRANSPOSE) {
144ff1e7120SSebastian Grimberg     // L-vector -> E-vector
145cf8cbdd6SSebastian Grimberg     CeedInt elem_size;
146cf8cbdd6SSebastian Grimberg 
147cf8cbdd6SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
148dce49693SSebastian Grimberg     const CeedInt block_size = elem_size < 1024 ? (elem_size > 32 ? elem_size : 32) : 1024;
149cf8cbdd6SSebastian Grimberg     const CeedInt grid       = CeedDivUpInt(impl->num_nodes, block_size);
15058549094SSebastian Grimberg 
151dce49693SSebastian Grimberg     switch (rstr_type) {
152dce49693SSebastian Grimberg       case CEED_RESTRICTION_STRIDED: {
153cf8cbdd6SSebastian Grimberg         void *args[] = {&d_u, &d_v};
15458549094SSebastian Grimberg 
155cf8cbdd6SSebastian Grimberg         CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args));
156dce49693SSebastian Grimberg       } break;
157b20a4af9SJeremy L Thompson       case CEED_RESTRICTION_POINTS:
158dce49693SSebastian Grimberg       case CEED_RESTRICTION_STANDARD: {
159a267acd1SJeremy L Thompson         void *args[] = {&impl->d_offsets, &d_u, &d_v};
160dce49693SSebastian Grimberg 
161cf8cbdd6SSebastian Grimberg         CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args));
162dce49693SSebastian Grimberg       } break;
163dce49693SSebastian Grimberg       case CEED_RESTRICTION_ORIENTED: {
164dce49693SSebastian Grimberg         if (use_signs) {
165a267acd1SJeremy L Thompson           void *args[] = {&impl->d_offsets, &impl->d_orients, &d_u, &d_v};
166dce49693SSebastian Grimberg 
167cf8cbdd6SSebastian Grimberg           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args));
168dce49693SSebastian Grimberg         } else {
169a267acd1SJeremy L Thompson           void *args[] = {&impl->d_offsets, &d_u, &d_v};
170dce49693SSebastian Grimberg 
171cf8cbdd6SSebastian Grimberg           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedNoTranspose, grid, block_size, args));
172dce49693SSebastian Grimberg         }
173dce49693SSebastian Grimberg       } break;
174dce49693SSebastian Grimberg       case CEED_RESTRICTION_CURL_ORIENTED: {
175dce49693SSebastian Grimberg         if (use_signs && use_orients) {
176a267acd1SJeremy L Thompson           void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v};
177dce49693SSebastian Grimberg 
178cf8cbdd6SSebastian Grimberg           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args));
179dce49693SSebastian Grimberg         } else if (use_orients) {
180a267acd1SJeremy L Thompson           void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v};
181dce49693SSebastian Grimberg 
182cf8cbdd6SSebastian Grimberg           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedNoTranspose, grid, block_size, args));
183dce49693SSebastian Grimberg         } else {
184a267acd1SJeremy L Thompson           void *args[] = {&impl->d_offsets, &d_u, &d_v};
185dce49693SSebastian Grimberg 
186cf8cbdd6SSebastian Grimberg           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedNoTranspose, grid, block_size, args));
187dce49693SSebastian Grimberg         }
188dce49693SSebastian Grimberg       } break;
189ff1e7120SSebastian Grimberg     }
190ff1e7120SSebastian Grimberg   } else {
191ff1e7120SSebastian Grimberg     // E-vector -> L-vector
192cf8cbdd6SSebastian Grimberg     const bool    is_deterministic = impl->d_l_vec_indices != NULL;
193dce49693SSebastian Grimberg     const CeedInt block_size       = 32;
194cf8cbdd6SSebastian Grimberg     const CeedInt grid             = CeedDivUpInt(impl->num_nodes, block_size);
195ca735530SJeremy L Thompson 
196dce49693SSebastian Grimberg     switch (rstr_type) {
197dce49693SSebastian Grimberg       case CEED_RESTRICTION_STRIDED: {
198cf8cbdd6SSebastian Grimberg         void *args[] = {&d_u, &d_v};
199dce49693SSebastian Grimberg 
200cf8cbdd6SSebastian Grimberg         CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
201dce49693SSebastian Grimberg       } break;
2020b63de31SJeremy L Thompson       case CEED_RESTRICTION_POINTS: {
2030b63de31SJeremy L Thompson         if (!is_deterministic) {
2040b63de31SJeremy L Thompson           void *args[] = {&impl->d_offsets, &impl->d_points_per_elem, &d_u, &d_v};
2050b63de31SJeremy L Thompson 
2060b63de31SJeremy L Thompson           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
2070b63de31SJeremy L Thompson         } else {
2080b63de31SJeremy L Thompson           void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_points_per_elem, &impl->d_t_offsets, &d_u, &d_v};
2090b63de31SJeremy L Thompson 
2100b63de31SJeremy L Thompson           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
2110b63de31SJeremy L Thompson         }
2120b63de31SJeremy L Thompson       } break;
213dce49693SSebastian Grimberg       case CEED_RESTRICTION_STANDARD: {
214cf8cbdd6SSebastian Grimberg         if (!is_deterministic) {
215a267acd1SJeremy L Thompson           void *args[] = {&impl->d_offsets, &d_u, &d_v};
21658549094SSebastian Grimberg 
217cf8cbdd6SSebastian Grimberg           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
218ff1e7120SSebastian Grimberg         } else {
21958549094SSebastian Grimberg           void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
22058549094SSebastian Grimberg 
221cf8cbdd6SSebastian Grimberg           CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
22258549094SSebastian Grimberg         }
223dce49693SSebastian Grimberg       } break;
224dce49693SSebastian Grimberg       case CEED_RESTRICTION_ORIENTED: {
225dce49693SSebastian Grimberg         if (use_signs) {
226cf8cbdd6SSebastian Grimberg           if (!is_deterministic) {
227a267acd1SJeremy L Thompson             void *args[] = {&impl->d_offsets, &impl->d_orients, &d_u, &d_v};
22858549094SSebastian Grimberg 
229cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
230dce49693SSebastian Grimberg           } else {
2317aa91133SSebastian Grimberg             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_orients, &d_u, &d_v};
2327aa91133SSebastian Grimberg 
233cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
2347aa91133SSebastian Grimberg           }
2357aa91133SSebastian Grimberg         } else {
236cf8cbdd6SSebastian Grimberg           if (!is_deterministic) {
237a267acd1SJeremy L Thompson             void *args[] = {&impl->d_offsets, &d_u, &d_v};
238dce49693SSebastian Grimberg 
239cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args));
240dce49693SSebastian Grimberg           } else {
241dce49693SSebastian Grimberg             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
242dce49693SSebastian Grimberg 
243cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args));
244dce49693SSebastian Grimberg           }
245dce49693SSebastian Grimberg         }
246dce49693SSebastian Grimberg       } break;
247dce49693SSebastian Grimberg       case CEED_RESTRICTION_CURL_ORIENTED: {
248dce49693SSebastian Grimberg         if (use_signs && use_orients) {
249cf8cbdd6SSebastian Grimberg           if (!is_deterministic) {
250a267acd1SJeremy L Thompson             void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v};
251dce49693SSebastian Grimberg 
252cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
2537aa91133SSebastian Grimberg           } else {
2547aa91133SSebastian Grimberg             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v};
2557aa91133SSebastian Grimberg 
256cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args));
2577aa91133SSebastian Grimberg           }
258dce49693SSebastian Grimberg         } else if (use_orients) {
259cf8cbdd6SSebastian Grimberg           if (!is_deterministic) {
260a267acd1SJeremy L Thompson             void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v};
261dce49693SSebastian Grimberg 
262cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args));
263dce49693SSebastian Grimberg           } else {
2647aa91133SSebastian Grimberg             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v};
2657aa91133SSebastian Grimberg 
266cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args));
2677aa91133SSebastian Grimberg           }
2687aa91133SSebastian Grimberg         } else {
269cf8cbdd6SSebastian Grimberg           if (!is_deterministic) {
270a267acd1SJeremy L Thompson             void *args[] = {&impl->d_offsets, &d_u, &d_v};
271dce49693SSebastian Grimberg 
272cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedTranspose, grid, block_size, args));
273dce49693SSebastian Grimberg           } else {
274dce49693SSebastian Grimberg             void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
275dce49693SSebastian Grimberg 
276cf8cbdd6SSebastian Grimberg             CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedTranspose, grid, block_size, args));
277dce49693SSebastian Grimberg           }
278dce49693SSebastian Grimberg         }
279dce49693SSebastian Grimberg       } break;
280ff1e7120SSebastian Grimberg     }
281ff1e7120SSebastian Grimberg   }
282ff1e7120SSebastian Grimberg 
283ff1e7120SSebastian Grimberg   if (request != CEED_REQUEST_IMMEDIATE && request != CEED_REQUEST_ORDERED) *request = NULL;
284ff1e7120SSebastian Grimberg 
285ff1e7120SSebastian Grimberg   // Restore arrays
286ff1e7120SSebastian Grimberg   CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
287ff1e7120SSebastian Grimberg   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
288ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
289ff1e7120SSebastian Grimberg }
290ff1e7120SSebastian Grimberg 
291ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
292dce49693SSebastian Grimberg // Apply restriction
293dce49693SSebastian Grimberg //------------------------------------------------------------------------------
294dce49693SSebastian Grimberg static int CeedElemRestrictionApply_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, CeedRequest *request) {
295dce49693SSebastian Grimberg   return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, true, true, u, v, request);
296dce49693SSebastian Grimberg }
297dce49693SSebastian Grimberg 
298dce49693SSebastian Grimberg //------------------------------------------------------------------------------
299dce49693SSebastian Grimberg // Apply unsigned restriction
300dce49693SSebastian Grimberg //------------------------------------------------------------------------------
301dce49693SSebastian Grimberg static int CeedElemRestrictionApplyUnsigned_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v,
302dce49693SSebastian Grimberg                                                  CeedRequest *request) {
303dce49693SSebastian Grimberg   return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, false, true, u, v, request);
304dce49693SSebastian Grimberg }
305dce49693SSebastian Grimberg 
306dce49693SSebastian Grimberg //------------------------------------------------------------------------------
307dce49693SSebastian Grimberg // Apply unoriented restriction
308dce49693SSebastian Grimberg //------------------------------------------------------------------------------
309dce49693SSebastian Grimberg static int CeedElemRestrictionApplyUnoriented_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v,
310dce49693SSebastian Grimberg                                                    CeedRequest *request) {
311dce49693SSebastian Grimberg   return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, false, false, u, v, request);
312dce49693SSebastian Grimberg }
313dce49693SSebastian Grimberg 
314dce49693SSebastian Grimberg //------------------------------------------------------------------------------
315ff1e7120SSebastian Grimberg // Get offsets
316ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
317ff1e7120SSebastian Grimberg static int CeedElemRestrictionGetOffsets_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt **offsets) {
318ff1e7120SSebastian Grimberg   CeedElemRestriction_Cuda *impl;
319b20a4af9SJeremy L Thompson   CeedRestrictionType       rstr_type;
320ff1e7120SSebastian Grimberg 
321ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
322b20a4af9SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type));
323ff1e7120SSebastian Grimberg   switch (mem_type) {
324ff1e7120SSebastian Grimberg     case CEED_MEM_HOST:
325b20a4af9SJeremy L Thompson       *offsets = rstr_type == CEED_RESTRICTION_POINTS ? impl->h_offsets_at_points : impl->h_offsets;
326ff1e7120SSebastian Grimberg       break;
327ff1e7120SSebastian Grimberg     case CEED_MEM_DEVICE:
328b20a4af9SJeremy L Thompson       *offsets = rstr_type == CEED_RESTRICTION_POINTS ? impl->d_offsets_at_points : impl->d_offsets;
329ff1e7120SSebastian Grimberg       break;
330ff1e7120SSebastian Grimberg   }
331ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
332ff1e7120SSebastian Grimberg }
333ff1e7120SSebastian Grimberg 
334ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
335dce49693SSebastian Grimberg // Get orientations
336dce49693SSebastian Grimberg //------------------------------------------------------------------------------
337dce49693SSebastian Grimberg static int CeedElemRestrictionGetOrientations_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const bool **orients) {
338dce49693SSebastian Grimberg   CeedElemRestriction_Cuda *impl;
339dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
340dce49693SSebastian Grimberg 
341dce49693SSebastian Grimberg   switch (mem_type) {
342dce49693SSebastian Grimberg     case CEED_MEM_HOST:
343dce49693SSebastian Grimberg       *orients = impl->h_orients;
344dce49693SSebastian Grimberg       break;
345dce49693SSebastian Grimberg     case CEED_MEM_DEVICE:
346dce49693SSebastian Grimberg       *orients = impl->d_orients;
347dce49693SSebastian Grimberg       break;
348dce49693SSebastian Grimberg   }
349dce49693SSebastian Grimberg   return CEED_ERROR_SUCCESS;
350dce49693SSebastian Grimberg }
351dce49693SSebastian Grimberg 
352dce49693SSebastian Grimberg //------------------------------------------------------------------------------
353dce49693SSebastian Grimberg // Get curl-conforming orientations
354dce49693SSebastian Grimberg //------------------------------------------------------------------------------
355dce49693SSebastian Grimberg static int CeedElemRestrictionGetCurlOrientations_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt8 **curl_orients) {
356dce49693SSebastian Grimberg   CeedElemRestriction_Cuda *impl;
357dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
358dce49693SSebastian Grimberg 
359dce49693SSebastian Grimberg   switch (mem_type) {
360dce49693SSebastian Grimberg     case CEED_MEM_HOST:
361dce49693SSebastian Grimberg       *curl_orients = impl->h_curl_orients;
362dce49693SSebastian Grimberg       break;
363dce49693SSebastian Grimberg     case CEED_MEM_DEVICE:
364dce49693SSebastian Grimberg       *curl_orients = impl->d_curl_orients;
365dce49693SSebastian Grimberg       break;
366dce49693SSebastian Grimberg   }
367dce49693SSebastian Grimberg   return CEED_ERROR_SUCCESS;
368dce49693SSebastian Grimberg }
369dce49693SSebastian Grimberg 
370dce49693SSebastian Grimberg //------------------------------------------------------------------------------
371b20a4af9SJeremy L Thompson // Get offset for padded AtPoints E-layout
372b20a4af9SJeremy L Thompson //------------------------------------------------------------------------------
373b20a4af9SJeremy L Thompson static int CeedElemRestrictionGetAtPointsElementOffset_Cuda(CeedElemRestriction rstr, CeedInt elem, CeedSize *elem_offset) {
374b20a4af9SJeremy L Thompson   CeedInt layout[3];
375b20a4af9SJeremy L Thompson 
376b20a4af9SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetELayout(rstr, layout));
377b20a4af9SJeremy L Thompson   *elem_offset = 0 * layout[0] + 0 * layout[1] + elem * layout[2];
378b20a4af9SJeremy L Thompson   return CEED_ERROR_SUCCESS;
379b20a4af9SJeremy L Thompson }
380b20a4af9SJeremy L Thompson 
381b20a4af9SJeremy L Thompson //------------------------------------------------------------------------------
382ff1e7120SSebastian Grimberg // Destroy restriction
383ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
384dce49693SSebastian Grimberg static int CeedElemRestrictionDestroy_Cuda(CeedElemRestriction rstr) {
385ff1e7120SSebastian Grimberg   Ceed                      ceed;
386ca735530SJeremy L Thompson   CeedElemRestriction_Cuda *impl;
387ca735530SJeremy L Thompson 
388dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
389dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
390cf8cbdd6SSebastian Grimberg   if (impl->module) {
391ff1e7120SSebastian Grimberg     CeedCallCuda(ceed, cuModuleUnload(impl->module));
392cf8cbdd6SSebastian Grimberg   }
393a267acd1SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->h_offsets_owned));
394f5d1e504SJeremy L Thompson   CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_offsets_owned));
395081aa29dSJeremy L Thompson   CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_t_offsets));
396081aa29dSJeremy L Thompson   CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_t_indices));
397081aa29dSJeremy L Thompson   CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_l_vec_indices));
398a267acd1SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->h_orients_owned));
399f5d1e504SJeremy L Thompson   CeedCallCuda(ceed, cudaFree((bool *)impl->d_orients_owned));
400a267acd1SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->h_curl_orients_owned));
401f5d1e504SJeremy L Thompson   CeedCallCuda(ceed, cudaFree((CeedInt8 *)impl->d_curl_orients_owned));
402b20a4af9SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->h_offsets_at_points_owned));
403b20a4af9SJeremy L Thompson   CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_offsets_at_points_owned));
4040b63de31SJeremy L Thompson   CeedCallBackend(CeedFree(&impl->h_points_per_elem_owned));
4050b63de31SJeremy L Thompson   CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_points_per_elem_owned));
406ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&impl));
407ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
408ff1e7120SSebastian Grimberg }
409ff1e7120SSebastian Grimberg 
410ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
411ff1e7120SSebastian Grimberg // Create transpose offsets and indices
412ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
413b20a4af9SJeremy L Thompson static int CeedElemRestrictionOffset_Cuda(const CeedElemRestriction rstr, const CeedInt elem_size, const CeedInt *indices) {
414ff1e7120SSebastian Grimberg   Ceed                      ceed;
415ca735530SJeremy L Thompson   bool                     *is_node;
416ff1e7120SSebastian Grimberg   CeedSize                  l_size;
417b20a4af9SJeremy L Thompson   CeedInt                   num_elem, num_comp, num_nodes = 0;
418ca735530SJeremy L Thompson   CeedInt                  *ind_to_offset, *l_vec_indices, *t_offsets, *t_indices;
419b20a4af9SJeremy L Thompson   CeedRestrictionType       rstr_type;
420ca735530SJeremy L Thompson   CeedElemRestriction_Cuda *impl;
421ca735530SJeremy L Thompson 
422dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
423dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
424dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem));
425b20a4af9SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type));
426dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size));
427dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp));
428ca735530SJeremy L Thompson   const CeedInt size_indices = num_elem * elem_size;
429ff1e7120SSebastian Grimberg 
430ff1e7120SSebastian Grimberg   // Count num_nodes
431ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(l_size, &is_node));
432ca735530SJeremy L Thompson 
433ff1e7120SSebastian Grimberg   for (CeedInt i = 0; i < size_indices; i++) is_node[indices[i]] = 1;
434ff1e7120SSebastian Grimberg   for (CeedInt i = 0; i < l_size; i++) num_nodes += is_node[i];
435ff1e7120SSebastian Grimberg   impl->num_nodes = num_nodes;
436ff1e7120SSebastian Grimberg 
437ff1e7120SSebastian Grimberg   // L-vector offsets array
438ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(l_size, &ind_to_offset));
439ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(num_nodes, &l_vec_indices));
440ca735530SJeremy L Thompson   for (CeedInt i = 0, j = 0; i < l_size; i++) {
441ff1e7120SSebastian Grimberg     if (is_node[i]) {
442ff1e7120SSebastian Grimberg       l_vec_indices[j] = i;
443ff1e7120SSebastian Grimberg       ind_to_offset[i] = j++;
444ff1e7120SSebastian Grimberg     }
445ff1e7120SSebastian Grimberg   }
446ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&is_node));
447ff1e7120SSebastian Grimberg 
448ff1e7120SSebastian Grimberg   // Compute transpose offsets and indices
449ff1e7120SSebastian Grimberg   const CeedInt size_offsets = num_nodes + 1;
450ca735530SJeremy L Thompson 
451ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(size_offsets, &t_offsets));
452ff1e7120SSebastian Grimberg   CeedCallBackend(CeedMalloc(size_indices, &t_indices));
453ff1e7120SSebastian Grimberg   // Count node multiplicity
454ff1e7120SSebastian Grimberg   for (CeedInt e = 0; e < num_elem; ++e) {
455ff1e7120SSebastian Grimberg     for (CeedInt i = 0; i < elem_size; ++i) ++t_offsets[ind_to_offset[indices[elem_size * e + i]] + 1];
456ff1e7120SSebastian Grimberg   }
457ff1e7120SSebastian Grimberg   // Convert to running sum
458ff1e7120SSebastian Grimberg   for (CeedInt i = 1; i < size_offsets; ++i) t_offsets[i] += t_offsets[i - 1];
459ff1e7120SSebastian Grimberg   // List all E-vec indices associated with L-vec node
460ff1e7120SSebastian Grimberg   for (CeedInt e = 0; e < num_elem; ++e) {
461ff1e7120SSebastian Grimberg     for (CeedInt i = 0; i < elem_size; ++i) {
462ff1e7120SSebastian Grimberg       const CeedInt lid = elem_size * e + i;
463ff1e7120SSebastian Grimberg       const CeedInt gid = indices[lid];
464ca735530SJeremy L Thompson 
465ff1e7120SSebastian Grimberg       t_indices[t_offsets[ind_to_offset[gid]]++] = lid;
466ff1e7120SSebastian Grimberg     }
467ff1e7120SSebastian Grimberg   }
468ff1e7120SSebastian Grimberg   // Reset running sum
469ff1e7120SSebastian Grimberg   for (int i = size_offsets - 1; i > 0; --i) t_offsets[i] = t_offsets[i - 1];
470ff1e7120SSebastian Grimberg   t_offsets[0] = 0;
471ff1e7120SSebastian Grimberg 
472ff1e7120SSebastian Grimberg   // Copy data to device
473ff1e7120SSebastian Grimberg   // -- L-vector indices
474ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_l_vec_indices, num_nodes * sizeof(CeedInt)));
475081aa29dSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_l_vec_indices, l_vec_indices, num_nodes * sizeof(CeedInt), cudaMemcpyHostToDevice));
476ff1e7120SSebastian Grimberg   // -- Transpose offsets
477ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_offsets, size_offsets * sizeof(CeedInt)));
478081aa29dSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_t_offsets, t_offsets, size_offsets * sizeof(CeedInt), cudaMemcpyHostToDevice));
479ff1e7120SSebastian Grimberg   // -- Transpose indices
480ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_indices, size_indices * sizeof(CeedInt)));
481081aa29dSJeremy L Thompson   CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_t_indices, t_indices, size_indices * sizeof(CeedInt), cudaMemcpyHostToDevice));
482ff1e7120SSebastian Grimberg 
483ff1e7120SSebastian Grimberg   // Cleanup
484ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&ind_to_offset));
485ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&l_vec_indices));
486ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&t_offsets));
487ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&t_indices));
488ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
489ff1e7120SSebastian Grimberg }
490ff1e7120SSebastian Grimberg 
491ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
492ff1e7120SSebastian Grimberg // Create restriction
493ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
494a267acd1SJeremy L Thompson int CeedElemRestrictionCreate_Cuda(CeedMemType mem_type, CeedCopyMode copy_mode, const CeedInt *offsets, const bool *orients,
495dce49693SSebastian Grimberg                                    const CeedInt8 *curl_orients, CeedElemRestriction rstr) {
496ca735530SJeremy L Thompson   Ceed                      ceed, ceed_parent;
497dce49693SSebastian Grimberg   bool                      is_deterministic;
498b20a4af9SJeremy L Thompson   CeedInt                   num_elem, num_comp, elem_size;
499ca735530SJeremy L Thompson   CeedRestrictionType       rstr_type;
500ff1e7120SSebastian Grimberg   CeedElemRestriction_Cuda *impl;
501ca735530SJeremy L Thompson 
502dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed));
503ca735530SJeremy L Thompson   CeedCallBackend(CeedGetParent(ceed, &ceed_parent));
504ca735530SJeremy L Thompson   CeedCallBackend(CeedIsDeterministic(ceed_parent, &is_deterministic));
505dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem));
506b20a4af9SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp));
507dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size));
50822eb1385SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type));
509b20a4af9SJeremy L Thompson   // Use max number of points as elem size for AtPoints restrictions
510b20a4af9SJeremy L Thompson   if (rstr_type == CEED_RESTRICTION_POINTS) {
511b20a4af9SJeremy L Thompson     CeedInt max_points = 0;
512b20a4af9SJeremy L Thompson 
513b20a4af9SJeremy L Thompson     for (CeedInt i = 0; i < num_elem; i++) {
514b20a4af9SJeremy L Thompson       max_points = CeedIntMax(max_points, offsets[i + 1] - offsets[i]);
515b20a4af9SJeremy L Thompson     }
516b20a4af9SJeremy L Thompson     elem_size = max_points;
517b20a4af9SJeremy L Thompson   }
518ca735530SJeremy L Thompson   const CeedInt size = num_elem * elem_size;
519ff1e7120SSebastian Grimberg 
520dce49693SSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &impl));
521dce49693SSebastian Grimberg   impl->num_nodes = size;
522dce49693SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionSetData(rstr, impl));
52322eb1385SJeremy L Thompson 
52422eb1385SJeremy L Thompson   // Set layouts
52522eb1385SJeremy L Thompson   {
52622eb1385SJeremy L Thompson     bool    has_backend_strides;
52722eb1385SJeremy L Thompson     CeedInt layout[3] = {1, size, elem_size};
52822eb1385SJeremy L Thompson 
529dce49693SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionSetELayout(rstr, layout));
53022eb1385SJeremy L Thompson     if (rstr_type == CEED_RESTRICTION_STRIDED) {
53122eb1385SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &has_backend_strides));
53222eb1385SJeremy L Thompson       if (has_backend_strides) {
53322eb1385SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionSetLLayout(rstr, layout));
53422eb1385SJeremy L Thompson       }
53522eb1385SJeremy L Thompson     }
53622eb1385SJeremy L Thompson   }
537ff1e7120SSebastian Grimberg 
538b20a4af9SJeremy L Thompson   // Pad AtPoints indices
539b20a4af9SJeremy L Thompson   if (rstr_type == CEED_RESTRICTION_POINTS) {
540b20a4af9SJeremy L Thompson     CeedSize offsets_len = elem_size * num_elem, at_points_size = num_elem + 1;
5410b63de31SJeremy L Thompson     CeedInt  max_points = elem_size, *offsets_padded, *points_per_elem;
542b20a4af9SJeremy L Thompson 
543b20a4af9SJeremy L Thompson     CeedCheck(mem_type == CEED_MEM_HOST, ceed, CEED_ERROR_BACKEND, "only MemType Host supported when creating AtPoints restriction");
544b20a4af9SJeremy L Thompson     CeedCallBackend(CeedMalloc(offsets_len, &offsets_padded));
5450b63de31SJeremy L Thompson     CeedCallBackend(CeedMalloc(num_elem, &points_per_elem));
546b20a4af9SJeremy L Thompson     for (CeedInt i = 0; i < num_elem; i++) {
547b20a4af9SJeremy L Thompson       CeedInt num_points = offsets[i + 1] - offsets[i];
548b20a4af9SJeremy L Thompson 
5490b63de31SJeremy L Thompson       points_per_elem[i] = num_points;
550b20a4af9SJeremy L Thompson       at_points_size += num_points;
551b20a4af9SJeremy L Thompson       // -- Copy all points in element
552b20a4af9SJeremy L Thompson       for (CeedInt j = 0; j < num_points; j++) {
5538be297eeSJeremy L Thompson         offsets_padded[i * max_points + j] = offsets[offsets[i] + j] * num_comp;
554b20a4af9SJeremy L Thompson       }
555b20a4af9SJeremy L Thompson       // -- Replicate out last point in element
556b20a4af9SJeremy L Thompson       for (CeedInt j = num_points; j < max_points; j++) {
5578be297eeSJeremy L Thompson         offsets_padded[i * max_points + j] = offsets[offsets[i] + num_points - 1] * num_comp;
558b20a4af9SJeremy L Thompson       }
559b20a4af9SJeremy L Thompson     }
560b20a4af9SJeremy L Thompson     CeedCallBackend(CeedSetHostCeedIntArray(offsets, copy_mode, at_points_size, &impl->h_offsets_at_points_owned, &impl->h_offsets_at_points_borrowed,
561b20a4af9SJeremy L Thompson                                             &impl->h_offsets_at_points));
562b20a4af9SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_offsets_at_points_owned, at_points_size * sizeof(CeedInt)));
563b20a4af9SJeremy L Thompson     CeedCallCuda(ceed, cudaMemcpy((CeedInt **)impl->d_offsets_at_points_owned, impl->h_offsets_at_points, at_points_size * sizeof(CeedInt),
564b20a4af9SJeremy L Thompson                                   cudaMemcpyHostToDevice));
565b20a4af9SJeremy L Thompson     impl->d_offsets_at_points = (CeedInt *)impl->d_offsets_at_points_owned;
566b20a4af9SJeremy L Thompson 
567b20a4af9SJeremy L Thompson     // -- Use padded offsets for the rest of the setup
568b20a4af9SJeremy L Thompson     offsets   = (const CeedInt *)offsets_padded;
569b20a4af9SJeremy L Thompson     copy_mode = CEED_OWN_POINTER;
5702e88d319SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionSetAtPointsEVectorSize(rstr, elem_size * num_elem * num_comp));
5710b63de31SJeremy L Thompson 
5720b63de31SJeremy L Thompson     // -- Points per element
5730b63de31SJeremy L Thompson     CeedCallBackend(CeedSetHostCeedIntArray(points_per_elem, CEED_OWN_POINTER, num_elem, &impl->h_points_per_elem_owned,
5740b63de31SJeremy L Thompson                                             &impl->h_points_per_elem_borrowed, &impl->h_points_per_elem));
5750b63de31SJeremy L Thompson     CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_points_per_elem_owned, num_elem * sizeof(CeedInt)));
5760b63de31SJeremy L Thompson     CeedCallCuda(ceed,
5770b63de31SJeremy L Thompson                  cudaMemcpy((CeedInt **)impl->d_points_per_elem_owned, impl->h_points_per_elem, num_elem * sizeof(CeedInt), cudaMemcpyHostToDevice));
5780b63de31SJeremy L Thompson     impl->d_points_per_elem = (CeedInt *)impl->d_points_per_elem_owned;
579b20a4af9SJeremy L Thompson   }
580b20a4af9SJeremy L Thompson 
581dce49693SSebastian Grimberg   // Set up device offset/orientation arrays
582dce49693SSebastian Grimberg   if (rstr_type != CEED_RESTRICTION_STRIDED) {
583ff1e7120SSebastian Grimberg     switch (mem_type) {
584ff1e7120SSebastian Grimberg       case CEED_MEM_HOST: {
585f5d1e504SJeremy L Thompson         CeedCallBackend(CeedSetHostCeedIntArray(offsets, copy_mode, size, &impl->h_offsets_owned, &impl->h_offsets_borrowed, &impl->h_offsets));
586a267acd1SJeremy L Thompson         CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_offsets_owned, size * sizeof(CeedInt)));
587f5d1e504SJeremy L Thompson         CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_offsets_owned, impl->h_offsets, size * sizeof(CeedInt), cudaMemcpyHostToDevice));
588f5d1e504SJeremy L Thompson         impl->d_offsets = (CeedInt *)impl->d_offsets_owned;
589b20a4af9SJeremy L Thompson         if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(rstr, elem_size, offsets));
590dce49693SSebastian Grimberg       } break;
591ff1e7120SSebastian Grimberg       case CEED_MEM_DEVICE: {
592f5d1e504SJeremy L Thompson         CeedCallBackend(CeedSetDeviceCeedIntArray_Cuda(ceed, offsets, copy_mode, size, &impl->d_offsets_owned, &impl->d_offsets_borrowed,
593f5d1e504SJeremy L Thompson                                                        (const CeedInt **)&impl->d_offsets));
594a267acd1SJeremy L Thompson         CeedCallBackend(CeedMalloc(size, &impl->h_offsets_owned));
595f5d1e504SJeremy L Thompson         CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->h_offsets_owned, impl->d_offsets, size * sizeof(CeedInt), cudaMemcpyDeviceToHost));
596a267acd1SJeremy L Thompson         impl->h_offsets = impl->h_offsets_owned;
597b20a4af9SJeremy L Thompson         if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(rstr, elem_size, offsets));
598dce49693SSebastian Grimberg       } break;
599ff1e7120SSebastian Grimberg     }
600ff1e7120SSebastian Grimberg 
601dce49693SSebastian Grimberg     // Orientation data
602dce49693SSebastian Grimberg     if (rstr_type == CEED_RESTRICTION_ORIENTED) {
603dce49693SSebastian Grimberg       switch (mem_type) {
604dce49693SSebastian Grimberg         case CEED_MEM_HOST: {
605f5d1e504SJeremy L Thompson           CeedCallBackend(CeedSetHostBoolArray(orients, copy_mode, size, &impl->h_orients_owned, &impl->h_orients_borrowed, &impl->h_orients));
606a267acd1SJeremy L Thompson           CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_orients_owned, size * sizeof(bool)));
607f5d1e504SJeremy L Thompson           CeedCallCuda(ceed, cudaMemcpy((bool *)impl->d_orients_owned, impl->h_orients, size * sizeof(bool), cudaMemcpyHostToDevice));
608a267acd1SJeremy L Thompson           impl->d_orients = impl->d_orients_owned;
609dce49693SSebastian Grimberg         } break;
610dce49693SSebastian Grimberg         case CEED_MEM_DEVICE: {
611f5d1e504SJeremy L Thompson           CeedCallBackend(CeedSetDeviceBoolArray_Cuda(ceed, orients, copy_mode, size, &impl->d_orients_owned, &impl->d_orients_borrowed,
612f5d1e504SJeremy L Thompson                                                       (const bool **)&impl->d_orients));
613a267acd1SJeremy L Thompson           CeedCallBackend(CeedMalloc(size, &impl->h_orients_owned));
614f5d1e504SJeremy L Thompson           CeedCallCuda(ceed, cudaMemcpy((bool *)impl->h_orients_owned, impl->d_orients, size * sizeof(bool), cudaMemcpyDeviceToHost));
615a267acd1SJeremy L Thompson           impl->h_orients = impl->h_orients_owned;
616dce49693SSebastian Grimberg         } break;
617dce49693SSebastian Grimberg       }
618dce49693SSebastian Grimberg     } else if (rstr_type == CEED_RESTRICTION_CURL_ORIENTED) {
619dce49693SSebastian Grimberg       switch (mem_type) {
620dce49693SSebastian Grimberg         case CEED_MEM_HOST: {
621f5d1e504SJeremy L Thompson           CeedCallBackend(CeedSetHostCeedInt8Array(curl_orients, copy_mode, 3 * size, &impl->h_curl_orients_owned, &impl->h_curl_orients_borrowed,
622f5d1e504SJeremy L Thompson                                                    &impl->h_curl_orients));
623a267acd1SJeremy L Thompson           CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_curl_orients_owned, 3 * size * sizeof(CeedInt8)));
624f5d1e504SJeremy L Thompson           CeedCallCuda(ceed,
625f5d1e504SJeremy L Thompson                        cudaMemcpy((CeedInt8 *)impl->d_curl_orients_owned, impl->h_curl_orients, 3 * size * sizeof(CeedInt8), cudaMemcpyHostToDevice));
626a267acd1SJeremy L Thompson           impl->d_curl_orients = impl->d_curl_orients_owned;
627dce49693SSebastian Grimberg         } break;
628dce49693SSebastian Grimberg         case CEED_MEM_DEVICE: {
629f5d1e504SJeremy L Thompson           CeedCallBackend(CeedSetDeviceCeedInt8Array_Cuda(ceed, curl_orients, copy_mode, 3 * size, &impl->d_curl_orients_owned,
630f5d1e504SJeremy L Thompson                                                           &impl->d_curl_orients_borrowed, (const CeedInt8 **)&impl->d_curl_orients));
631a267acd1SJeremy L Thompson           CeedCallBackend(CeedMalloc(3 * size, &impl->h_curl_orients_owned));
632f5d1e504SJeremy L Thompson           CeedCallCuda(ceed,
633f5d1e504SJeremy L Thompson                        cudaMemcpy((CeedInt8 *)impl->h_curl_orients_owned, impl->d_curl_orients, 3 * size * sizeof(CeedInt8), cudaMemcpyDeviceToHost));
634a267acd1SJeremy L Thompson           impl->h_curl_orients = impl->h_curl_orients_owned;
635dce49693SSebastian Grimberg         } break;
636dce49693SSebastian Grimberg       }
637dce49693SSebastian Grimberg     }
638dce49693SSebastian Grimberg   }
639ca735530SJeremy L Thompson 
640ff1e7120SSebastian Grimberg   // Register backend functions
641dce49693SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Apply", CeedElemRestrictionApply_Cuda));
642dce49693SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnsigned", CeedElemRestrictionApplyUnsigned_Cuda));
643dce49693SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnoriented", CeedElemRestrictionApplyUnoriented_Cuda));
644dce49693SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOffsets", CeedElemRestrictionGetOffsets_Cuda));
645dce49693SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOrientations", CeedElemRestrictionGetOrientations_Cuda));
646dce49693SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetCurlOrientations", CeedElemRestrictionGetCurlOrientations_Cuda));
647b20a4af9SJeremy L Thompson   if (rstr_type == CEED_RESTRICTION_POINTS) {
648b20a4af9SJeremy L Thompson     CeedCallBackend(
649b20a4af9SJeremy L Thompson         CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetAtPointsElementOffset", CeedElemRestrictionGetAtPointsElementOffset_Cuda));
650b20a4af9SJeremy L Thompson   }
651dce49693SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Destroy", CeedElemRestrictionDestroy_Cuda));
652ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
653ff1e7120SSebastian Grimberg }
654ff1e7120SSebastian Grimberg 
655ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
656