xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-ref/ceed-cuda-ref-restriction.c (revision ca7355308a14805110a7d98dabf1f658d5ec16d5)
1ff1e7120SSebastian Grimberg // Copyright (c) 2017-2022, 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 //------------------------------------------------------------------------------
22ff1e7120SSebastian Grimberg // Apply restriction
23ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
24ff1e7120SSebastian Grimberg static int CeedElemRestrictionApply_Cuda(CeedElemRestriction r, CeedTransposeMode t_mode, CeedVector u, CeedVector v, CeedRequest *request) {
25ff1e7120SSebastian Grimberg   Ceed                      ceed;
26ff1e7120SSebastian Grimberg   Ceed_Cuda                *data;
27ff1e7120SSebastian Grimberg   CUfunction                kernel;
28*ca735530SJeremy L Thompson   CeedInt                   num_elem, elem_size;
29ff1e7120SSebastian Grimberg   const CeedScalar         *d_u;
30ff1e7120SSebastian Grimberg   CeedScalar               *d_v;
31*ca735530SJeremy L Thompson   CeedElemRestriction_Cuda *impl;
32*ca735530SJeremy L Thompson 
33*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(r, &impl));
34*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
35*ca735530SJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &data));
36*ca735530SJeremy L Thompson   CeedElemRestrictionGetNumElements(r, &num_elem);
37*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size));
38*ca735530SJeremy L Thompson   const CeedInt num_nodes = impl->num_nodes;
39*ca735530SJeremy L Thompson 
40*ca735530SJeremy L Thompson   // Get vectors
41ff1e7120SSebastian Grimberg   CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
42ff1e7120SSebastian Grimberg   if (t_mode == CEED_TRANSPOSE) {
43ff1e7120SSebastian Grimberg     // Sum into for transpose mode, e-vec to l-vec
44ff1e7120SSebastian Grimberg     CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v));
45ff1e7120SSebastian Grimberg   } else {
46ff1e7120SSebastian Grimberg     // Overwrite for notranspose mode, l-vec to e-vec
47ff1e7120SSebastian Grimberg     CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
48ff1e7120SSebastian Grimberg   }
49ff1e7120SSebastian Grimberg 
50ff1e7120SSebastian Grimberg   // Restrict
51ff1e7120SSebastian Grimberg   if (t_mode == CEED_NOTRANSPOSE) {
52ff1e7120SSebastian Grimberg     // L-vector -> E-vector
53ff1e7120SSebastian Grimberg     if (impl->d_ind) {
54ff1e7120SSebastian Grimberg       // -- Offsets provided
55ff1e7120SSebastian Grimberg       kernel             = impl->OffsetNoTranspose;
56ff1e7120SSebastian Grimberg       void   *args[]     = {&num_elem, &impl->d_ind, &d_u, &d_v};
57ff1e7120SSebastian Grimberg       CeedInt block_size = elem_size < 1024 ? (elem_size > 32 ? elem_size : 32) : 1024;
5858549094SSebastian Grimberg 
59ff1e7120SSebastian Grimberg       CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
60ff1e7120SSebastian Grimberg     } else {
61ff1e7120SSebastian Grimberg       // -- Strided restriction
62ff1e7120SSebastian Grimberg       kernel             = impl->StridedNoTranspose;
63ff1e7120SSebastian Grimberg       void   *args[]     = {&num_elem, &d_u, &d_v};
64ff1e7120SSebastian Grimberg       CeedInt block_size = elem_size < 1024 ? (elem_size > 32 ? elem_size : 32) : 1024;
6558549094SSebastian Grimberg 
66ff1e7120SSebastian Grimberg       CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
67ff1e7120SSebastian Grimberg     }
68ff1e7120SSebastian Grimberg   } else {
69ff1e7120SSebastian Grimberg     // E-vector -> L-vector
70ff1e7120SSebastian Grimberg     if (impl->d_ind) {
71ff1e7120SSebastian Grimberg       // -- Offsets provided
7258549094SSebastian Grimberg       CeedInt block_size = 32;
73*ca735530SJeremy L Thompson 
7458549094SSebastian Grimberg       if (impl->OffsetTranspose) {
75ff1e7120SSebastian Grimberg         kernel       = impl->OffsetTranspose;
7658549094SSebastian Grimberg         void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v};
7758549094SSebastian Grimberg 
78ff1e7120SSebastian Grimberg         CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
79ff1e7120SSebastian Grimberg       } else {
8058549094SSebastian Grimberg         kernel       = impl->OffsetTransposeDet;
8158549094SSebastian Grimberg         void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v};
8258549094SSebastian Grimberg 
8358549094SSebastian Grimberg         CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
8458549094SSebastian Grimberg       }
8558549094SSebastian Grimberg     } else {
86ff1e7120SSebastian Grimberg       // -- Strided restriction
87ff1e7120SSebastian Grimberg       kernel             = impl->StridedTranspose;
88ff1e7120SSebastian Grimberg       void   *args[]     = {&num_elem, &d_u, &d_v};
8958549094SSebastian Grimberg       CeedInt block_size = 32;
9058549094SSebastian Grimberg 
91ff1e7120SSebastian Grimberg       CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args));
92ff1e7120SSebastian Grimberg     }
93ff1e7120SSebastian Grimberg   }
94ff1e7120SSebastian Grimberg 
95ff1e7120SSebastian Grimberg   if (request != CEED_REQUEST_IMMEDIATE && request != CEED_REQUEST_ORDERED) *request = NULL;
96ff1e7120SSebastian Grimberg 
97ff1e7120SSebastian Grimberg   // Restore arrays
98ff1e7120SSebastian Grimberg   CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
99ff1e7120SSebastian Grimberg   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
100ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
101ff1e7120SSebastian Grimberg }
102ff1e7120SSebastian Grimberg 
103ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
104ff1e7120SSebastian Grimberg // Get offsets
105ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
106ff1e7120SSebastian Grimberg static int CeedElemRestrictionGetOffsets_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt **offsets) {
107ff1e7120SSebastian Grimberg   CeedElemRestriction_Cuda *impl;
108ff1e7120SSebastian Grimberg 
109*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl));
110ff1e7120SSebastian Grimberg   switch (mem_type) {
111ff1e7120SSebastian Grimberg     case CEED_MEM_HOST:
112ff1e7120SSebastian Grimberg       *offsets = impl->h_ind;
113ff1e7120SSebastian Grimberg       break;
114ff1e7120SSebastian Grimberg     case CEED_MEM_DEVICE:
115ff1e7120SSebastian Grimberg       *offsets = impl->d_ind;
116ff1e7120SSebastian Grimberg       break;
117ff1e7120SSebastian Grimberg   }
118ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
119ff1e7120SSebastian Grimberg }
120ff1e7120SSebastian Grimberg 
121ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
122ff1e7120SSebastian Grimberg // Destroy restriction
123ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
124ff1e7120SSebastian Grimberg static int CeedElemRestrictionDestroy_Cuda(CeedElemRestriction r) {
125ff1e7120SSebastian Grimberg   Ceed                      ceed;
126*ca735530SJeremy L Thompson   CeedElemRestriction_Cuda *impl;
127*ca735530SJeremy L Thompson 
128*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(r, &impl));
129ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
130ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cuModuleUnload(impl->module));
131ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&impl->h_ind_allocated));
132ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaFree(impl->d_ind_allocated));
133ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaFree(impl->d_t_offsets));
134ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaFree(impl->d_t_indices));
135ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaFree(impl->d_l_vec_indices));
136ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&impl));
137ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
138ff1e7120SSebastian Grimberg }
139ff1e7120SSebastian Grimberg 
140ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
141ff1e7120SSebastian Grimberg // Create transpose offsets and indices
142ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
143ff1e7120SSebastian Grimberg static int CeedElemRestrictionOffset_Cuda(const CeedElemRestriction r, const CeedInt *indices) {
144ff1e7120SSebastian Grimberg   Ceed                      ceed;
145*ca735530SJeremy L Thompson   bool                     *is_node;
146ff1e7120SSebastian Grimberg   CeedSize                  l_size;
147*ca735530SJeremy L Thompson   CeedInt                   num_elem, elem_size, num_comp, num_nodes = 0;
148*ca735530SJeremy L Thompson   CeedInt                  *ind_to_offset, *l_vec_indices, *t_offsets, *t_indices;
149*ca735530SJeremy L Thompson   CeedElemRestriction_Cuda *impl;
150*ca735530SJeremy L Thompson 
151*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
152*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetData(r, &impl));
153ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(r, &num_elem));
154ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size));
155ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetLVectorSize(r, &l_size));
156ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(r, &num_comp));
157*ca735530SJeremy L Thompson   const CeedInt size_indices = num_elem * elem_size;
158ff1e7120SSebastian Grimberg 
159ff1e7120SSebastian Grimberg   // Count num_nodes
160ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(l_size, &is_node));
161*ca735530SJeremy L Thompson 
162ff1e7120SSebastian Grimberg   for (CeedInt i = 0; i < size_indices; i++) is_node[indices[i]] = 1;
163ff1e7120SSebastian Grimberg   for (CeedInt i = 0; i < l_size; i++) num_nodes += is_node[i];
164ff1e7120SSebastian Grimberg   impl->num_nodes = num_nodes;
165ff1e7120SSebastian Grimberg 
166ff1e7120SSebastian Grimberg   // L-vector offsets array
167ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(l_size, &ind_to_offset));
168ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(num_nodes, &l_vec_indices));
169*ca735530SJeremy L Thompson   for (CeedInt i = 0, j = 0; i < l_size; i++) {
170ff1e7120SSebastian Grimberg     if (is_node[i]) {
171ff1e7120SSebastian Grimberg       l_vec_indices[j] = i;
172ff1e7120SSebastian Grimberg       ind_to_offset[i] = j++;
173ff1e7120SSebastian Grimberg     }
174ff1e7120SSebastian Grimberg   }
175ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&is_node));
176ff1e7120SSebastian Grimberg 
177ff1e7120SSebastian Grimberg   // Compute transpose offsets and indices
178ff1e7120SSebastian Grimberg   const CeedInt size_offsets = num_nodes + 1;
179*ca735530SJeremy L Thompson 
180ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(size_offsets, &t_offsets));
181ff1e7120SSebastian Grimberg   CeedCallBackend(CeedMalloc(size_indices, &t_indices));
182ff1e7120SSebastian Grimberg   // Count node multiplicity
183ff1e7120SSebastian Grimberg   for (CeedInt e = 0; e < num_elem; ++e) {
184ff1e7120SSebastian Grimberg     for (CeedInt i = 0; i < elem_size; ++i) ++t_offsets[ind_to_offset[indices[elem_size * e + i]] + 1];
185ff1e7120SSebastian Grimberg   }
186ff1e7120SSebastian Grimberg   // Convert to running sum
187ff1e7120SSebastian Grimberg   for (CeedInt i = 1; i < size_offsets; ++i) t_offsets[i] += t_offsets[i - 1];
188ff1e7120SSebastian Grimberg   // List all E-vec indices associated with L-vec node
189ff1e7120SSebastian Grimberg   for (CeedInt e = 0; e < num_elem; ++e) {
190ff1e7120SSebastian Grimberg     for (CeedInt i = 0; i < elem_size; ++i) {
191ff1e7120SSebastian Grimberg       const CeedInt lid = elem_size * e + i;
192ff1e7120SSebastian Grimberg       const CeedInt gid = indices[lid];
193*ca735530SJeremy L Thompson 
194ff1e7120SSebastian Grimberg       t_indices[t_offsets[ind_to_offset[gid]]++] = lid;
195ff1e7120SSebastian Grimberg     }
196ff1e7120SSebastian Grimberg   }
197ff1e7120SSebastian Grimberg   // Reset running sum
198ff1e7120SSebastian Grimberg   for (int i = size_offsets - 1; i > 0; --i) t_offsets[i] = t_offsets[i - 1];
199ff1e7120SSebastian Grimberg   t_offsets[0] = 0;
200ff1e7120SSebastian Grimberg 
201ff1e7120SSebastian Grimberg   // Copy data to device
202ff1e7120SSebastian Grimberg   // -- L-vector indices
203ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_l_vec_indices, num_nodes * sizeof(CeedInt)));
204ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(impl->d_l_vec_indices, l_vec_indices, num_nodes * sizeof(CeedInt), cudaMemcpyHostToDevice));
205ff1e7120SSebastian Grimberg   // -- Transpose offsets
206ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_offsets, size_offsets * sizeof(CeedInt)));
207ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(impl->d_t_offsets, t_offsets, size_offsets * sizeof(CeedInt), cudaMemcpyHostToDevice));
208ff1e7120SSebastian Grimberg   // -- Transpose indices
209ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_indices, size_indices * sizeof(CeedInt)));
210ff1e7120SSebastian Grimberg   CeedCallCuda(ceed, cudaMemcpy(impl->d_t_indices, t_indices, size_indices * sizeof(CeedInt), cudaMemcpyHostToDevice));
211ff1e7120SSebastian Grimberg 
212ff1e7120SSebastian Grimberg   // Cleanup
213ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&ind_to_offset));
214ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&l_vec_indices));
215ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&t_offsets));
216ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&t_indices));
217ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
218ff1e7120SSebastian Grimberg }
219ff1e7120SSebastian Grimberg 
220ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
221ff1e7120SSebastian Grimberg // Create restriction
222ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
223fcbe8c06SSebastian Grimberg int CeedElemRestrictionCreate_Cuda(CeedMemType mem_type, CeedCopyMode copy_mode, const CeedInt *indices, const bool *orients,
2240c73c039SSebastian Grimberg                                    const CeedInt8 *curl_orients, CeedElemRestriction r) {
225*ca735530SJeremy L Thompson   Ceed                      ceed, ceed_parent;
226*ca735530SJeremy L Thompson   bool                      is_deterministic, is_strided;
227*ca735530SJeremy L Thompson   CeedInt                   num_elem, num_comp, elem_size, comp_stride = 1;
228*ca735530SJeremy L Thompson   CeedRestrictionType       rstr_type;
229ff1e7120SSebastian Grimberg   CeedElemRestriction_Cuda *impl;
230*ca735530SJeremy L Thompson 
231*ca735530SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed));
232ff1e7120SSebastian Grimberg   CeedCallBackend(CeedCalloc(1, &impl));
233*ca735530SJeremy L Thompson   CeedCallBackend(CeedGetParent(ceed, &ceed_parent));
234*ca735530SJeremy L Thompson   CeedCallBackend(CeedIsDeterministic(ceed_parent, &is_deterministic));
235ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumElements(r, &num_elem));
236ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetNumComponents(r, &num_comp));
237ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size));
238*ca735530SJeremy L Thompson   const CeedInt size       = num_elem * elem_size;
239ff1e7120SSebastian Grimberg   CeedInt       strides[3] = {1, size, elem_size};
240*ca735530SJeremy L Thompson   CeedInt       layout[3]  = {1, elem_size * num_elem, elem_size};
241ff1e7120SSebastian Grimberg 
24200125730SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionGetType(r, &rstr_type));
24300125730SSebastian Grimberg   CeedCheck(rstr_type != CEED_RESTRICTION_ORIENTED && rstr_type != CEED_RESTRICTION_CURL_ORIENTED, ceed, CEED_ERROR_BACKEND,
24400125730SSebastian Grimberg             "Backend does not implement CeedElemRestrictionCreateOriented or CeedElemRestrictionCreateCurlOriented");
24500125730SSebastian Grimberg 
246ff1e7120SSebastian Grimberg   // Stride data
247ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionIsStrided(r, &is_strided));
248ff1e7120SSebastian Grimberg   if (is_strided) {
249ff1e7120SSebastian Grimberg     bool has_backend_strides;
250*ca735530SJeremy L Thompson 
251ff1e7120SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionHasBackendStrides(r, &has_backend_strides));
252ff1e7120SSebastian Grimberg     if (!has_backend_strides) {
253ff1e7120SSebastian Grimberg       CeedCallBackend(CeedElemRestrictionGetStrides(r, &strides));
254ff1e7120SSebastian Grimberg     }
255ff1e7120SSebastian Grimberg   } else {
256ff1e7120SSebastian Grimberg     CeedCallBackend(CeedElemRestrictionGetCompStride(r, &comp_stride));
257ff1e7120SSebastian Grimberg   }
258ff1e7120SSebastian Grimberg 
259ff1e7120SSebastian Grimberg   impl->h_ind           = NULL;
260ff1e7120SSebastian Grimberg   impl->h_ind_allocated = NULL;
261ff1e7120SSebastian Grimberg   impl->d_ind           = NULL;
262ff1e7120SSebastian Grimberg   impl->d_ind_allocated = NULL;
263ff1e7120SSebastian Grimberg   impl->d_t_indices     = NULL;
264ff1e7120SSebastian Grimberg   impl->d_t_offsets     = NULL;
265ff1e7120SSebastian Grimberg   impl->num_nodes       = size;
266ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionSetData(r, impl));
267ff1e7120SSebastian Grimberg   CeedCallBackend(CeedElemRestrictionSetELayout(r, layout));
268ff1e7120SSebastian Grimberg 
269ff1e7120SSebastian Grimberg   // Set up device indices/offset arrays
270ff1e7120SSebastian Grimberg   switch (mem_type) {
271ff1e7120SSebastian Grimberg     case CEED_MEM_HOST: {
272ff1e7120SSebastian Grimberg       switch (copy_mode) {
273ff1e7120SSebastian Grimberg         case CEED_OWN_POINTER:
274ff1e7120SSebastian Grimberg           impl->h_ind_allocated = (CeedInt *)indices;
275ff1e7120SSebastian Grimberg           impl->h_ind           = (CeedInt *)indices;
276ff1e7120SSebastian Grimberg           break;
277ff1e7120SSebastian Grimberg         case CEED_USE_POINTER:
278ff1e7120SSebastian Grimberg           impl->h_ind = (CeedInt *)indices;
279ff1e7120SSebastian Grimberg           break;
280ff1e7120SSebastian Grimberg         case CEED_COPY_VALUES:
281ff1e7120SSebastian Grimberg           if (indices != NULL) {
282ff1e7120SSebastian Grimberg             CeedCallBackend(CeedMalloc(elem_size * num_elem, &impl->h_ind_allocated));
283ff1e7120SSebastian Grimberg             memcpy(impl->h_ind_allocated, indices, elem_size * num_elem * sizeof(CeedInt));
284ff1e7120SSebastian Grimberg             impl->h_ind = impl->h_ind_allocated;
285ff1e7120SSebastian Grimberg           }
286ff1e7120SSebastian Grimberg           break;
287ff1e7120SSebastian Grimberg       }
288ff1e7120SSebastian Grimberg       if (indices != NULL) {
289ff1e7120SSebastian Grimberg         CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_ind, size * sizeof(CeedInt)));
290ff1e7120SSebastian Grimberg         impl->d_ind_allocated = impl->d_ind;  // We own the device memory
291ff1e7120SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), cudaMemcpyHostToDevice));
29258549094SSebastian Grimberg         if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(r, indices));
293ff1e7120SSebastian Grimberg       }
294ff1e7120SSebastian Grimberg       break;
295ff1e7120SSebastian Grimberg     }
296ff1e7120SSebastian Grimberg     case CEED_MEM_DEVICE: {
297ff1e7120SSebastian Grimberg       switch (copy_mode) {
298ff1e7120SSebastian Grimberg         case CEED_COPY_VALUES:
299ff1e7120SSebastian Grimberg           if (indices != NULL) {
300ff1e7120SSebastian Grimberg             CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_ind, size * sizeof(CeedInt)));
301ff1e7120SSebastian Grimberg             impl->d_ind_allocated = impl->d_ind;  // We own the device memory
302ff1e7120SSebastian Grimberg             CeedCallCuda(ceed, cudaMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), cudaMemcpyDeviceToDevice));
303ff1e7120SSebastian Grimberg           }
304ff1e7120SSebastian Grimberg           break;
305ff1e7120SSebastian Grimberg         case CEED_OWN_POINTER:
306ff1e7120SSebastian Grimberg           impl->d_ind           = (CeedInt *)indices;
307ff1e7120SSebastian Grimberg           impl->d_ind_allocated = impl->d_ind;
308ff1e7120SSebastian Grimberg           break;
309ff1e7120SSebastian Grimberg         case CEED_USE_POINTER:
310ff1e7120SSebastian Grimberg           impl->d_ind = (CeedInt *)indices;
311ff1e7120SSebastian Grimberg       }
312ff1e7120SSebastian Grimberg       if (indices != NULL) {
313ff1e7120SSebastian Grimberg         CeedCallBackend(CeedMalloc(elem_size * num_elem, &impl->h_ind_allocated));
314ff1e7120SSebastian Grimberg         CeedCallCuda(ceed, cudaMemcpy(impl->h_ind_allocated, impl->d_ind, elem_size * num_elem * sizeof(CeedInt), cudaMemcpyDeviceToHost));
315ff1e7120SSebastian Grimberg         impl->h_ind = impl->h_ind_allocated;
31658549094SSebastian Grimberg         if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(r, indices));
317ff1e7120SSebastian Grimberg       }
318ff1e7120SSebastian Grimberg       break;
319ff1e7120SSebastian Grimberg     }
320ff1e7120SSebastian Grimberg     // LCOV_EXCL_START
321ff1e7120SSebastian Grimberg     default:
322ff1e7120SSebastian Grimberg       return CeedError(ceed, CEED_ERROR_BACKEND, "Only MemType = HOST or DEVICE supported");
323ff1e7120SSebastian Grimberg       // LCOV_EXCL_STOP
324ff1e7120SSebastian Grimberg   }
325ff1e7120SSebastian Grimberg 
32658549094SSebastian Grimberg   // Compile CUDA kernels (add atomicAdd function for old NVidia architectures)
327ff1e7120SSebastian Grimberg   CeedInt num_nodes = impl->num_nodes;
32858549094SSebastian Grimberg   char   *restriction_kernel_path, *restriction_kernel_source = NULL;
329*ca735530SJeremy L Thompson 
330ff1e7120SSebastian Grimberg   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction.h", &restriction_kernel_path));
33123d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n");
33258549094SSebastian Grimberg   if (!is_deterministic) {
33358549094SSebastian Grimberg     struct cudaDeviceProp prop;
33458549094SSebastian Grimberg     Ceed_Cuda            *ceed_data;
335*ca735530SJeremy L Thompson 
33658549094SSebastian Grimberg     CeedCallBackend(CeedGetData(ceed, &ceed_data));
33758549094SSebastian Grimberg     CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id));
33858549094SSebastian Grimberg     if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) {
33958549094SSebastian Grimberg       char *atomic_add_path;
340*ca735530SJeremy L Thompson 
34158549094SSebastian Grimberg       CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-atomic-add-fallback.h", &atomic_add_path));
34258549094SSebastian Grimberg       CeedCallBackend(CeedLoadSourceToBuffer(ceed, atomic_add_path, &restriction_kernel_source));
34358549094SSebastian Grimberg       CeedCallBackend(CeedLoadSourceToInitializedBuffer(ceed, restriction_kernel_path, &restriction_kernel_source));
34458549094SSebastian Grimberg       CeedCallBackend(CeedFree(&atomic_add_path));
34558549094SSebastian Grimberg     }
34658549094SSebastian Grimberg   }
34758549094SSebastian Grimberg   if (!restriction_kernel_source) {
348ff1e7120SSebastian Grimberg     CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source));
34958549094SSebastian Grimberg   }
35023d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n");
351*ca735530SJeremy L Thompson   CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 8, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem,
352*ca735530SJeremy L Thompson                                    "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", num_nodes, "RSTR_COMP_STRIDE", comp_stride, "RSTR_STRIDE_NODES",
353*ca735530SJeremy L Thompson                                    strides[0], "RSTR_STRIDE_COMP", strides[1], "RSTR_STRIDE_ELEM", strides[2]));
354ff1e7120SSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedNoTranspose", &impl->StridedNoTranspose));
35558549094SSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedTranspose", &impl->StridedTranspose));
356ff1e7120SSebastian Grimberg   CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->OffsetNoTranspose));
35758549094SSebastian Grimberg   if (!is_deterministic) {
35858549094SSebastian Grimberg     CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->OffsetTranspose));
35958549094SSebastian Grimberg   } else {
36058549094SSebastian Grimberg     CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTransposeDet", &impl->OffsetTransposeDet));
36158549094SSebastian Grimberg   }
362ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&restriction_kernel_path));
363ff1e7120SSebastian Grimberg   CeedCallBackend(CeedFree(&restriction_kernel_source));
364ff1e7120SSebastian Grimberg 
365ff1e7120SSebastian Grimberg   // Register backend functions
366ff1e7120SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "Apply", CeedElemRestrictionApply_Cuda));
367ff1e7120SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "ApplyUnsigned", CeedElemRestrictionApply_Cuda));
3687c1dbaffSSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "ApplyUnoriented", CeedElemRestrictionApply_Cuda));
369ff1e7120SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "GetOffsets", CeedElemRestrictionGetOffsets_Cuda));
370ff1e7120SSebastian Grimberg   CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "Destroy", CeedElemRestrictionDestroy_Cuda));
371ff1e7120SSebastian Grimberg   return CEED_ERROR_SUCCESS;
372ff1e7120SSebastian Grimberg }
373ff1e7120SSebastian Grimberg 
374ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------
375